xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
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 
23d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
24d71ae5a4SJacob 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);
353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36af33a6ddSJed Brown }
37af33a6ddSJed Brown 
38d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicDestroy(Characteristic *c)
39d71ae5a4SJacob Faibussowitsch {
40af33a6ddSJed Brown   PetscFunctionBegin;
413ba16761SJacob Faibussowitsch   if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
426bf464f9SBarry Smith   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
433ba16761SJacob Faibussowitsch   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
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));
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59af33a6ddSJed Brown }
60af33a6ddSJed Brown 
61d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
62d71ae5a4SJacob 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;
1083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
109af33a6ddSJed Brown }
110af33a6ddSJed Brown 
111af33a6ddSJed Brown /*@C
112af33a6ddSJed Brown    CharacteristicSetType - Builds Characteristic for a particular solver.
113af33a6ddSJed Brown 
114c3339decSBarry Smith    Logically Collective
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 
124bcf0153eSBarry Smith   Level: intermediate
125bcf0153eSBarry Smith 
126af33a6ddSJed Brown    Notes:
12730a76a96SBarry Smith    See "include/petsccharacteristic.h" for available methods
128af33a6ddSJed Brown 
129af33a6ddSJed Brown   Normally, it is best to use the CharacteristicSetFromOptions() command and
130af33a6ddSJed Brown   then set the Characteristic type from the options database rather than by using
131af33a6ddSJed Brown   this routine.  Using the options database provides the user with
132af33a6ddSJed Brown   maximum flexibility in evaluating the many different Krylov methods.
133af33a6ddSJed Brown   The CharacteristicSetType() routine is provided for those situations where it
134af33a6ddSJed Brown   is necessary to set the iterative solver independently of the command
135af33a6ddSJed Brown   line or options database.  This might be the case, for example, when
136af33a6ddSJed Brown   the choice of iterative solver changes during the execution of the
137af33a6ddSJed Brown   program, and the user's application is taking responsibility for
138af33a6ddSJed Brown   choosing the appropriate method.  In other words, this routine is
139af33a6ddSJed Brown   not for beginners.
140af33a6ddSJed Brown 
141*1cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicType`
142af33a6ddSJed Brown @*/
143d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
144d71ae5a4SJacob Faibussowitsch {
145af33a6ddSJed Brown   PetscBool match;
1465f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(Characteristic);
147af33a6ddSJed Brown 
148af33a6ddSJed Brown   PetscFunctionBegin;
149af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
150af33a6ddSJed Brown   PetscValidCharPointer(type, 2);
151af33a6ddSJed Brown 
1529566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)c, type, &match));
1533ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
154af33a6ddSJed Brown 
155af33a6ddSJed Brown   if (c->data) {
156af33a6ddSJed Brown     /* destroy the old private Characteristic context */
157dbbe0bcdSBarry Smith     PetscUseTypeMethod(c, destroy);
1580298fd71SBarry Smith     c->ops->destroy = NULL;
1592120983fSLisandro Dalcin     c->data         = NULL;
160af33a6ddSJed Brown   }
161af33a6ddSJed Brown 
1629566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(CharacteristicList, type, &r));
1633c633725SBarry Smith   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
164af33a6ddSJed Brown   c->setupcalled = 0;
1659566063dSJacob Faibussowitsch   PetscCall((*r)(c));
1669566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)c, type));
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168af33a6ddSJed Brown }
169af33a6ddSJed Brown 
170af33a6ddSJed Brown /*@
171af33a6ddSJed Brown    CharacteristicSetUp - Sets up the internal data structures for the
172af33a6ddSJed Brown    later use of an iterative solver.
173af33a6ddSJed Brown 
174c3339decSBarry Smith    Collective
175af33a6ddSJed Brown 
176af33a6ddSJed Brown    Input Parameter:
177af33a6ddSJed Brown .  ksp   - iterative context obtained from CharacteristicCreate()
178af33a6ddSJed Brown 
179af33a6ddSJed Brown    Level: developer
180af33a6ddSJed Brown 
181*1cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicCreate()`, `CharacteristicSolve()`, `CharacteristicDestroy()`
182af33a6ddSJed Brown @*/
183d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetUp(Characteristic c)
184d71ae5a4SJacob Faibussowitsch {
185af33a6ddSJed Brown   PetscFunctionBegin;
186af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
187af33a6ddSJed Brown 
18848a46eb9SPierre Jolivet   if (!((PetscObject)c)->type_name) PetscCall(CharacteristicSetType(c, CHARACTERISTICDA));
189af33a6ddSJed Brown 
1903ba16761SJacob Faibussowitsch   if (c->setupcalled == 2) PetscFunctionReturn(PETSC_SUCCESS);
191af33a6ddSJed Brown 
1929566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
193ad540459SPierre Jolivet   if (!c->setupcalled) PetscUseTypeMethod(c, setup);
1949566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
195af33a6ddSJed Brown   c->setupcalled = 2;
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
197af33a6ddSJed Brown }
198af33a6ddSJed Brown 
199af33a6ddSJed Brown /*@C
2001c84c290SBarry Smith   CharacteristicRegister -  Adds a solver to the method of characteristics package.
2011c84c290SBarry Smith 
2021c84c290SBarry Smith    Not Collective
2031c84c290SBarry Smith 
2041c84c290SBarry Smith    Input Parameters:
2052fe279fdSBarry Smith +  sname - name of a new user-defined solver
2062fe279fdSBarry Smith -  function - routine to create method context
2071c84c290SBarry Smith 
208bcf0153eSBarry Smith   Level: advanced
209bcf0153eSBarry Smith 
2101c84c290SBarry Smith   Sample usage:
2111c84c290SBarry Smith .vb
212bdf89e91SBarry Smith     CharacteristicRegister("my_char", MyCharCreate);
2131c84c290SBarry Smith .ve
2141c84c290SBarry Smith 
2151c84c290SBarry Smith   Then, your Characteristic type can be chosen with the procedural interface via
2161c84c290SBarry Smith .vb
2171c84c290SBarry Smith     CharacteristicCreate(MPI_Comm, Characteristic* &char);
2181c84c290SBarry Smith     CharacteristicSetType(char,"my_char");
2191c84c290SBarry Smith .ve
2201c84c290SBarry Smith    or at runtime via the option
2211c84c290SBarry Smith .vb
2221c84c290SBarry Smith     -characteristic_type my_char
2231c84c290SBarry Smith .ve
2241c84c290SBarry Smith 
2251c84c290SBarry Smith    Notes:
2261c84c290SBarry Smith    CharacteristicRegister() may be called multiple times to add several user-defined solvers.
2271c84c290SBarry Smith 
228*1cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicRegisterAll()`, `CharacteristicRegisterDestroy()`
229af33a6ddSJed Brown @*/
230d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicRegister(const char sname[], PetscErrorCode (*function)(Characteristic))
231d71ae5a4SJacob Faibussowitsch {
232af33a6ddSJed Brown   PetscFunctionBegin;
2339566063dSJacob Faibussowitsch   PetscCall(CharacteristicInitializePackage());
2349566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&CharacteristicList, sname, function));
2353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
236af33a6ddSJed Brown }
237af33a6ddSJed Brown 
238d71ae5a4SJacob 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)
239d71ae5a4SJacob Faibussowitsch {
240af33a6ddSJed Brown   PetscFunctionBegin;
241af33a6ddSJed Brown   c->velocityDA      = da;
242af33a6ddSJed Brown   c->velocity        = v;
243af33a6ddSJed Brown   c->velocityOld     = vOld;
244af33a6ddSJed Brown   c->numVelocityComp = numComponents;
245af33a6ddSJed Brown   c->velocityComp    = components;
246af33a6ddSJed Brown   c->velocityInterp  = interp;
247af33a6ddSJed Brown   c->velocityCtx     = ctx;
2483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
249af33a6ddSJed Brown }
250af33a6ddSJed Brown 
251d71ae5a4SJacob 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)
252d71ae5a4SJacob Faibussowitsch {
253af33a6ddSJed Brown   PetscFunctionBegin;
254af33a6ddSJed Brown   c->velocityDA          = da;
255af33a6ddSJed Brown   c->velocity            = v;
256af33a6ddSJed Brown   c->velocityOld         = vOld;
257af33a6ddSJed Brown   c->numVelocityComp     = numComponents;
258af33a6ddSJed Brown   c->velocityComp        = components;
259af33a6ddSJed Brown   c->velocityInterpLocal = interp;
260af33a6ddSJed Brown   c->velocityCtx         = ctx;
2613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
262af33a6ddSJed Brown }
263af33a6ddSJed Brown 
264d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx)
265d71ae5a4SJacob Faibussowitsch {
266af33a6ddSJed Brown   PetscFunctionBegin;
267af33a6ddSJed Brown #if 0
2683c633725SBarry 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.");
269af33a6ddSJed Brown #endif
270af33a6ddSJed Brown   c->fieldDA      = da;
271af33a6ddSJed Brown   c->field        = v;
272af33a6ddSJed Brown   c->numFieldComp = numComponents;
273af33a6ddSJed Brown   c->fieldComp    = components;
274af33a6ddSJed Brown   c->fieldInterp  = interp;
275af33a6ddSJed Brown   c->fieldCtx     = ctx;
2763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
277af33a6ddSJed Brown }
278af33a6ddSJed Brown 
279d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx)
280d71ae5a4SJacob Faibussowitsch {
281af33a6ddSJed Brown   PetscFunctionBegin;
282af33a6ddSJed Brown #if 0
2833c633725SBarry 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.");
284af33a6ddSJed Brown #endif
285af33a6ddSJed Brown   c->fieldDA          = da;
286af33a6ddSJed Brown   c->field            = v;
287af33a6ddSJed Brown   c->numFieldComp     = numComponents;
288af33a6ddSJed Brown   c->fieldComp        = components;
289af33a6ddSJed Brown   c->fieldInterpLocal = interp;
290af33a6ddSJed Brown   c->fieldCtx         = ctx;
2913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
292af33a6ddSJed Brown }
293af33a6ddSJed Brown 
294d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
295d71ae5a4SJacob Faibussowitsch {
296af33a6ddSJed Brown   CharacteristicPointDA2D Qi;
297af33a6ddSJed Brown   DM                      da = c->velocityDA;
298af33a6ddSJed Brown   Vec                     velocityLocal, velocityLocalOld;
299af33a6ddSJed Brown   Vec                     fieldLocal;
300af33a6ddSJed Brown   DMDALocalInfo           info;
301af33a6ddSJed Brown   PetscScalar           **solArray;
302af33a6ddSJed Brown   void                   *velocityArray;
303af33a6ddSJed Brown   void                   *velocityArrayOld;
304af33a6ddSJed Brown   void                   *fieldArray;
3051cc8b206SBarry Smith   PetscScalar            *interpIndices;
3061cc8b206SBarry Smith   PetscScalar            *velocityValues, *velocityValuesOld;
3071cc8b206SBarry Smith   PetscScalar            *fieldValues;
308af33a6ddSJed Brown   PetscMPIInt             rank;
309af33a6ddSJed Brown   PetscInt                dim;
310af33a6ddSJed Brown   PetscMPIInt             neighbors[9];
311af33a6ddSJed Brown   PetscInt                dof;
312af33a6ddSJed Brown   PetscInt                gx, gy;
313af33a6ddSJed Brown   PetscInt                n, is, ie, js, je, comp;
314af33a6ddSJed Brown 
315af33a6ddSJed Brown   PetscFunctionBegin;
316af33a6ddSJed Brown   c->queueSize = 0;
3179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
3189566063dSJacob Faibussowitsch   PetscCall(DMDAGetNeighborsRank(da, neighbors));
3199566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetNeighbors(c, 9, neighbors));
3209566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetUp(c));
321af33a6ddSJed Brown   /* global and local grid info */
3229566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, &gx, &gy, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3239566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
3249371c9d4SSatish Balay   is = info.xs;
3259371c9d4SSatish Balay   ie = info.xs + info.xm;
3269371c9d4SSatish Balay   js = info.ys;
3279371c9d4SSatish Balay   je = info.ys + info.ym;
328af33a6ddSJed Brown   /* Allocation */
3299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim, &interpIndices));
3309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValues));
3319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValuesOld));
3329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numFieldComp, &fieldValues));
3339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
334af33a6ddSJed Brown 
335af33a6ddSJed Brown   /* -----------------------------------------------------------------------
336af33a6ddSJed Brown      PART 1, AT t-dt/2
337af33a6ddSJed Brown      -----------------------------------------------------------------------*/
3389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
339af33a6ddSJed Brown   /* GET POSITION AT HALF TIME IN THE PAST */
340af33a6ddSJed Brown   if (c->velocityInterpLocal) {
3419566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocal));
3429566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocalOld));
3439566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3449566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3459566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3469566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3479566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray));
3489566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
349af33a6ddSJed Brown   }
3509566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n"));
351af33a6ddSJed Brown   for (Qi.j = js; Qi.j < je; Qi.j++) {
352af33a6ddSJed Brown     for (Qi.i = is; Qi.i < ie; Qi.i++) {
353af33a6ddSJed Brown       interpIndices[0] = Qi.i;
354af33a6ddSJed Brown       interpIndices[1] = Qi.j;
3559566063dSJacob Faibussowitsch       if (c->velocityInterpLocal) PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3569566063dSJacob Faibussowitsch       else PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
357af33a6ddSJed Brown       Qi.x = Qi.i - velocityValues[0] * dt / 2.0;
358af33a6ddSJed Brown       Qi.y = Qi.j - velocityValues[1] * dt / 2.0;
359af33a6ddSJed Brown 
360af33a6ddSJed Brown       /* Determine whether the position at t - dt/2 is local */
361af33a6ddSJed Brown       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
362af33a6ddSJed Brown 
363af33a6ddSJed Brown       /* Check for Periodic boundaries and move all periodic points back onto the domain */
3649566063dSJacob Faibussowitsch       PetscCall(DMDAMapCoordsToPeriodicDomain(da, &(Qi.x), &(Qi.y)));
3659566063dSJacob Faibussowitsch       PetscCall(CharacteristicAddPoint(c, &Qi));
366af33a6ddSJed Brown     }
367af33a6ddSJed Brown   }
3689566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
369af33a6ddSJed Brown 
3709566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
3719566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
3729566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
373af33a6ddSJed Brown 
3749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
375af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (local values) */
3769566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n"));
377af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
378af33a6ddSJed Brown     Qi = c->queue[n];
379af33a6ddSJed Brown     if (c->neighbors[Qi.proc] == rank) {
380af33a6ddSJed Brown       interpIndices[0] = Qi.x;
381af33a6ddSJed Brown       interpIndices[1] = Qi.y;
382af33a6ddSJed Brown       if (c->velocityInterpLocal) {
3839566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3849566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
385af33a6ddSJed Brown       } else {
3869566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3879566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
388af33a6ddSJed Brown       }
389af33a6ddSJed Brown       Qi.x = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
390af33a6ddSJed Brown       Qi.y = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
391af33a6ddSJed Brown     }
392af33a6ddSJed Brown     c->queue[n] = Qi;
393af33a6ddSJed Brown   }
3949566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
395af33a6ddSJed Brown 
3969566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
3979566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
3989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
399af33a6ddSJed Brown 
400af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (fill remote requests) */
4019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
40263a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote velocities at t_{n - 1/2}\n", c->queueRemoteSize));
403af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
404af33a6ddSJed Brown     Qi               = c->queueRemote[n];
405af33a6ddSJed Brown     interpIndices[0] = Qi.x;
406af33a6ddSJed Brown     interpIndices[1] = Qi.y;
407af33a6ddSJed Brown     if (c->velocityInterpLocal) {
4089566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4099566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
410af33a6ddSJed Brown     } else {
4119566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4129566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
413af33a6ddSJed Brown     }
414af33a6ddSJed Brown     Qi.x              = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
415af33a6ddSJed Brown     Qi.y              = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
416af33a6ddSJed Brown     c->queueRemote[n] = Qi;
417af33a6ddSJed Brown   }
4189566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
4199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
4209566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
4219566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
422af33a6ddSJed Brown   if (c->velocityInterpLocal) {
4239566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray));
4249566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
4259566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocal));
4269566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocalOld));
427af33a6ddSJed Brown   }
4289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
429af33a6ddSJed Brown 
430af33a6ddSJed Brown   /* -----------------------------------------------------------------------
431af33a6ddSJed Brown      PART 2, AT t-dt
432af33a6ddSJed Brown      -----------------------------------------------------------------------*/
433af33a6ddSJed Brown 
434af33a6ddSJed Brown   /* GET POSITION AT t_n (local values) */
4359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
4369566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n}\n"));
437af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
438af33a6ddSJed Brown     Qi   = c->queue[n];
439af33a6ddSJed Brown     Qi.x = Qi.i - Qi.x * dt;
440af33a6ddSJed Brown     Qi.y = Qi.j - Qi.y * dt;
441af33a6ddSJed Brown 
442af33a6ddSJed Brown     /* Determine whether the position at t-dt is local */
443af33a6ddSJed Brown     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
444af33a6ddSJed Brown 
445af33a6ddSJed Brown     /* Check for Periodic boundaries and move all periodic points back onto the domain */
4469566063dSJacob Faibussowitsch     PetscCall(DMDAMapCoordsToPeriodicDomain(da, &(Qi.x), &(Qi.y)));
447af33a6ddSJed Brown 
448af33a6ddSJed Brown     c->queue[n] = Qi;
449af33a6ddSJed Brown   }
4509566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
451af33a6ddSJed Brown 
4529566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4539566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
4549566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
455af33a6ddSJed Brown 
456af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
4579566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
458af33a6ddSJed Brown   if (c->fieldInterpLocal) {
4599566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->fieldDA, &fieldLocal));
4609566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4619566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4629566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray));
463af33a6ddSJed Brown   }
4649566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local field at t_{n}\n"));
465af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
466af33a6ddSJed Brown     if (c->neighbors[c->queue[n].proc] == rank) {
467af33a6ddSJed Brown       interpIndices[0] = c->queue[n].x;
468af33a6ddSJed Brown       interpIndices[1] = c->queue[n].y;
4699566063dSJacob Faibussowitsch       if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4709566063dSJacob Faibussowitsch       else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
471bbd56ea5SKarl Rupp       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
472af33a6ddSJed Brown     }
473af33a6ddSJed Brown   }
4749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
475af33a6ddSJed Brown 
4769566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4779566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
4789566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
479af33a6ddSJed Brown 
480af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
4819566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
48263a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote field points at t_{n}\n", c->queueRemoteSize));
483af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
484af33a6ddSJed Brown     interpIndices[0] = c->queueRemote[n].x;
485af33a6ddSJed Brown     interpIndices[1] = c->queueRemote[n].y;
486af33a6ddSJed Brown 
487af33a6ddSJed Brown     /* for debugging purposes */
488af33a6ddSJed Brown     if (1) { /* hacked bounds test...let's do better */
4899371c9d4SSatish Balay       PetscScalar im = interpIndices[0];
4909371c9d4SSatish Balay       PetscScalar jm = interpIndices[1];
491af33a6ddSJed Brown 
49263a3b9bcSJacob 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));
493af33a6ddSJed Brown     }
494af33a6ddSJed Brown 
4959566063dSJacob Faibussowitsch     if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4969566063dSJacob Faibussowitsch     else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
497bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
498af33a6ddSJed Brown   }
4999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
500af33a6ddSJed Brown 
5019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
5029566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
5039566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
504af33a6ddSJed Brown   if (c->fieldInterpLocal) {
5059566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray));
5069566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->fieldDA, &fieldLocal));
507af33a6ddSJed Brown   }
5089566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
509af33a6ddSJed Brown 
510af33a6ddSJed Brown   /* Return field of characteristics at t_n-1 */
5119566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5129566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(c->fieldDA, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5139566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(c->fieldDA, solution, &solArray));
514af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
515af33a6ddSJed Brown     Qi = c->queue[n];
516bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i * dof + c->fieldComp[comp]] = Qi.field[comp];
517af33a6ddSJed Brown   }
5189566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(c->fieldDA, solution, &solArray));
5199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5209566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
521af33a6ddSJed Brown 
522af33a6ddSJed Brown   /* Cleanup */
5239566063dSJacob Faibussowitsch   PetscCall(PetscFree(interpIndices));
5249566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValues));
5259566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValuesOld));
5269566063dSJacob Faibussowitsch   PetscCall(PetscFree(fieldValues));
5273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
528af33a6ddSJed Brown }
529af33a6ddSJed Brown 
530d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
531d71ae5a4SJacob Faibussowitsch {
532af33a6ddSJed Brown   PetscFunctionBegin;
533af33a6ddSJed Brown   c->numNeighbors = numNeighbors;
5349566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->neighbors));
5359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numNeighbors, &c->neighbors));
5369566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(c->neighbors, neighbors, numNeighbors));
5373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
538af33a6ddSJed Brown }
539af33a6ddSJed Brown 
540d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
541d71ae5a4SJacob Faibussowitsch {
542af33a6ddSJed Brown   PetscFunctionBegin;
54363a3b9bcSJacob Faibussowitsch   PetscCheck(c->queueSize < c->queueMax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exceeded maximum queue size %" PetscInt_FMT, c->queueMax);
544af33a6ddSJed Brown   c->queue[c->queueSize++] = *point;
5453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
546af33a6ddSJed Brown }
547af33a6ddSJed Brown 
5483ba16761SJacob Faibussowitsch PetscErrorCode CharacteristicSendCoordinatesBegin(Characteristic c)
549d71ae5a4SJacob Faibussowitsch {
550af33a6ddSJed Brown   PetscMPIInt rank, tag = 121;
551af33a6ddSJed Brown   PetscInt    i, n;
552af33a6ddSJed Brown 
553af33a6ddSJed Brown   PetscFunctionBegin;
5549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
5559566063dSJacob Faibussowitsch   PetscCall(CharacteristicHeapSort(c, c->queue, c->queueSize));
5569566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(c->needCount, c->numNeighbors));
557bbd56ea5SKarl Rupp   for (i = 0; i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
558af33a6ddSJed Brown   c->fillCount[0] = 0;
55948a46eb9SPierre 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])));
56048a46eb9SPierre Jolivet   for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
5619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
562af33a6ddSJed Brown   /* Initialize the remote queue */
563af33a6ddSJed Brown   c->queueLocalMax = c->localOffsets[0] = 0;
564af33a6ddSJed Brown   c->queueRemoteMax = c->remoteOffsets[0] = 0;
565af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
566af33a6ddSJed Brown     c->remoteOffsets[n] = c->queueRemoteMax;
567af33a6ddSJed Brown     c->queueRemoteMax += c->fillCount[n];
568af33a6ddSJed Brown     c->localOffsets[n] = c->queueLocalMax;
569af33a6ddSJed Brown     c->queueLocalMax += c->needCount[n];
570af33a6ddSJed Brown   }
571af33a6ddSJed Brown   /* HACK BEGIN */
572bbd56ea5SKarl Rupp   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
573af33a6ddSJed Brown   c->needCount[0] = 0;
574af33a6ddSJed Brown   /* HACK END */
575af33a6ddSJed Brown   if (c->queueRemoteMax) {
5769566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(c->queueRemoteMax, &c->queueRemote));
5770298fd71SBarry Smith   } else c->queueRemote = NULL;
578af33a6ddSJed Brown   c->queueRemoteSize = c->queueRemoteMax;
579af33a6ddSJed Brown 
580af33a6ddSJed Brown   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
581af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
58263a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Receiving %" PetscInt_FMT " requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]));
5839566063dSJacob 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])));
584af33a6ddSJed Brown   }
585af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
58663a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Sending %" PetscInt_FMT " requests for values from proc %d\n", c->needCount[n], c->neighbors[n]));
5879566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
588af33a6ddSJed Brown   }
5893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
590af33a6ddSJed Brown }
591af33a6ddSJed Brown 
592d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
593d71ae5a4SJacob Faibussowitsch {
594af33a6ddSJed Brown #if 0
595af33a6ddSJed Brown   PetscMPIInt rank;
596af33a6ddSJed Brown   PetscInt    n;
597af33a6ddSJed Brown #endif
598af33a6ddSJed Brown 
599af33a6ddSJed Brown   PetscFunctionBegin;
6009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
601af33a6ddSJed Brown #if 0
6029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
603af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
6043c633725SBarry 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);
605af33a6ddSJed Brown   }
606af33a6ddSJed Brown #endif
6073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
608af33a6ddSJed Brown }
609af33a6ddSJed Brown 
610d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
611d71ae5a4SJacob Faibussowitsch {
612af33a6ddSJed Brown   PetscMPIInt tag = 121;
613af33a6ddSJed Brown   PetscInt    n;
614af33a6ddSJed Brown 
615af33a6ddSJed Brown   PetscFunctionBegin;
616d5b43468SJose E. Roman   /* SEND AND RECEIVE FILLED REQUESTS for velocities at t_n+1/2 */
61748a46eb9SPierre 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])));
61848a46eb9SPierre 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)));
6193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
620af33a6ddSJed Brown }
621af33a6ddSJed Brown 
622d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
623d71ae5a4SJacob Faibussowitsch {
624af33a6ddSJed Brown   PetscFunctionBegin;
6259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
626af33a6ddSJed Brown   /* Free queue of requests from other procs */
6279566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->queueRemote));
6283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
629af33a6ddSJed Brown }
630af33a6ddSJed Brown 
631af33a6ddSJed Brown /*---------------------------------------------------------------------*/
632af33a6ddSJed Brown /*
633af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
634af33a6ddSJed Brown */
6350b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
636af33a6ddSJed Brown /*---------------------------------------------------------------------*/
637af33a6ddSJed Brown {
638af33a6ddSJed Brown   CharacteristicPointDA2D temp;
6390b8f38c8SBarry Smith   PetscInt                n;
640af33a6ddSJed Brown 
6410b8f38c8SBarry Smith   PetscFunctionBegin;
642af33a6ddSJed Brown   if (0) { /* Check the order of the queue before sorting */
6439566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Before Heap sort\n"));
64448a46eb9SPierre Jolivet     for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
645af33a6ddSJed Brown   }
646af33a6ddSJed Brown 
647af33a6ddSJed Brown   /* SORTING PHASE */
6489371c9d4SSatish Balay   for (n = (size / 2) - 1; n >= 0; n--) { PetscCall(CharacteristicSiftDown(c, queue, n, size - 1)); /* Rich had size-1 here, Matt had size*/ }
649af33a6ddSJed Brown   for (n = size - 1; n >= 1; n--) {
650af33a6ddSJed Brown     temp     = queue[0];
651af33a6ddSJed Brown     queue[0] = queue[n];
652af33a6ddSJed Brown     queue[n] = temp;
6539566063dSJacob Faibussowitsch     PetscCall(CharacteristicSiftDown(c, queue, 0, n - 1));
654af33a6ddSJed Brown   }
655af33a6ddSJed Brown   if (0) { /* Check the order of the queue after sorting */
6569566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Avter  Heap sort\n"));
65748a46eb9SPierre Jolivet     for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
658af33a6ddSJed Brown   }
6593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
660af33a6ddSJed Brown }
661af33a6ddSJed Brown 
662af33a6ddSJed Brown /*---------------------------------------------------------------------*/
663af33a6ddSJed Brown /*
664af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
665af33a6ddSJed Brown */
6660b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
667af33a6ddSJed Brown /*---------------------------------------------------------------------*/
668af33a6ddSJed Brown {
669af33a6ddSJed Brown   PetscBool               done = PETSC_FALSE;
670af33a6ddSJed Brown   PetscInt                maxChild;
671af33a6ddSJed Brown   CharacteristicPointDA2D temp;
672af33a6ddSJed Brown 
673af33a6ddSJed Brown   PetscFunctionBegin;
674af33a6ddSJed Brown   while ((root * 2 <= bottom) && (!done)) {
675af33a6ddSJed Brown     if (root * 2 == bottom) maxChild = root * 2;
676af33a6ddSJed Brown     else if (queue[root * 2].proc > queue[root * 2 + 1].proc) maxChild = root * 2;
677af33a6ddSJed Brown     else maxChild = root * 2 + 1;
678af33a6ddSJed Brown 
679af33a6ddSJed Brown     if (queue[root].proc < queue[maxChild].proc) {
680af33a6ddSJed Brown       temp            = queue[root];
681af33a6ddSJed Brown       queue[root]     = queue[maxChild];
682af33a6ddSJed Brown       queue[maxChild] = temp;
683af33a6ddSJed Brown       root            = maxChild;
684db4deed7SKarl Rupp     } else done = PETSC_TRUE;
685af33a6ddSJed Brown   }
6863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
687af33a6ddSJed Brown }
688af33a6ddSJed Brown 
689af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
690d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
691d71ae5a4SJacob Faibussowitsch {
692bff4a2f0SMatthew G. Knepley   DMBoundaryType bx, by;
693af33a6ddSJed Brown   PetscBool      IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
694af33a6ddSJed Brown   MPI_Comm       comm;
695af33a6ddSJed Brown   PetscMPIInt    rank;
696af33a6ddSJed Brown   PetscInt     **procs, pi, pj, pim, pip, pjm, pjp, PI, PJ;
697af33a6ddSJed Brown 
698af33a6ddSJed Brown   PetscFunctionBegin;
6999566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
7009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7019566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, &PI, &PJ, NULL, NULL, NULL, &bx, &by, NULL, NULL));
702af33a6ddSJed Brown 
703bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
704bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
705af33a6ddSJed Brown 
706af33a6ddSJed Brown   neighbors[0] = rank;
707af33a6ddSJed Brown   rank         = 0;
7089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(PJ, &procs));
709af33a6ddSJed Brown   for (pj = 0; pj < PJ; pj++) {
7109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PI, &(procs[pj])));
711af33a6ddSJed Brown     for (pi = 0; pi < PI; pi++) {
712af33a6ddSJed Brown       procs[pj][pi] = rank;
713af33a6ddSJed Brown       rank++;
714af33a6ddSJed Brown     }
715af33a6ddSJed Brown   }
716af33a6ddSJed Brown 
717af33a6ddSJed Brown   pi  = neighbors[0] % PI;
718af33a6ddSJed Brown   pj  = neighbors[0] / PI;
7199371c9d4SSatish Balay   pim = pi - 1;
7209371c9d4SSatish Balay   if (pim < 0) pim = PI - 1;
721af33a6ddSJed Brown   pip = (pi + 1) % PI;
7229371c9d4SSatish Balay   pjm = pj - 1;
7239371c9d4SSatish Balay   if (pjm < 0) pjm = PJ - 1;
724af33a6ddSJed Brown   pjp = (pj + 1) % PJ;
725af33a6ddSJed Brown 
726af33a6ddSJed Brown   neighbors[1] = procs[pj][pim];
727af33a6ddSJed Brown   neighbors[2] = procs[pjp][pim];
728af33a6ddSJed Brown   neighbors[3] = procs[pjp][pi];
729af33a6ddSJed Brown   neighbors[4] = procs[pjp][pip];
730af33a6ddSJed Brown   neighbors[5] = procs[pj][pip];
731af33a6ddSJed Brown   neighbors[6] = procs[pjm][pip];
732af33a6ddSJed Brown   neighbors[7] = procs[pjm][pi];
733af33a6ddSJed Brown   neighbors[8] = procs[pjm][pim];
734af33a6ddSJed Brown 
735af33a6ddSJed Brown   if (!IPeriodic) {
736af33a6ddSJed Brown     if (pi == 0) neighbors[1] = neighbors[2] = neighbors[8] = neighbors[0];
737af33a6ddSJed Brown     if (pi == PI - 1) neighbors[4] = neighbors[5] = neighbors[6] = neighbors[0];
738af33a6ddSJed Brown   }
739af33a6ddSJed Brown 
740af33a6ddSJed Brown   if (!JPeriodic) {
741af33a6ddSJed Brown     if (pj == 0) neighbors[6] = neighbors[7] = neighbors[8] = neighbors[0];
742af33a6ddSJed Brown     if (pj == PJ - 1) neighbors[2] = neighbors[3] = neighbors[4] = neighbors[0];
743af33a6ddSJed Brown   }
744af33a6ddSJed Brown 
74548a46eb9SPierre Jolivet   for (pj = 0; pj < PJ; pj++) PetscCall(PetscFree(procs[pj]));
7469566063dSJacob Faibussowitsch   PetscCall(PetscFree(procs));
7473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
748af33a6ddSJed Brown }
749af33a6ddSJed Brown 
750af33a6ddSJed Brown /*
751af33a6ddSJed Brown   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
752af33a6ddSJed Brown     2 | 3 | 4
753af33a6ddSJed Brown     __|___|__
754af33a6ddSJed Brown     1 | 0 | 5
755af33a6ddSJed Brown     __|___|__
756af33a6ddSJed Brown     8 | 7 | 6
757af33a6ddSJed Brown       |   |
758af33a6ddSJed Brown */
759d71ae5a4SJacob Faibussowitsch PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr)
760d71ae5a4SJacob Faibussowitsch {
761af33a6ddSJed Brown   DMDALocalInfo info;
7621cc8b206SBarry Smith   PetscReal     is, ie, js, je;
763af33a6ddSJed Brown 
7649566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
7659371c9d4SSatish Balay   is = (PetscReal)info.xs - 0.5;
7669371c9d4SSatish Balay   ie = (PetscReal)info.xs + info.xm - 0.5;
7679371c9d4SSatish Balay   js = (PetscReal)info.ys - 0.5;
7689371c9d4SSatish Balay   je = (PetscReal)info.ys + info.ym - 0.5;
769af33a6ddSJed Brown 
770af33a6ddSJed Brown   if (ir >= is && ir <= ie) { /* center column */
771bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 0;
772bbd56ea5SKarl Rupp     else if (jr < js) return 7;
773bbd56ea5SKarl Rupp     else return 3;
774af33a6ddSJed Brown   } else if (ir < is) { /* left column */
775bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 1;
776bbd56ea5SKarl Rupp     else if (jr < js) return 8;
777bbd56ea5SKarl Rupp     else return 2;
778af33a6ddSJed Brown   } else { /* right column */
779bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 5;
780bbd56ea5SKarl Rupp     else if (jr < js) return 6;
781bbd56ea5SKarl Rupp     else return 4;
782af33a6ddSJed Brown   }
783af33a6ddSJed Brown }
784