xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
239371c9d4SSatish Balay PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer) {
24af33a6ddSJed Brown   PetscBool iascii;
25af33a6ddSJed Brown 
26af33a6ddSJed Brown   PetscFunctionBegin;
27af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
28*48a46eb9SPierre 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));
339371c9d4SSatish Balay   if (!iascii) { PetscTryTypeMethod(c, view, viewer); }
34af33a6ddSJed Brown   PetscFunctionReturn(0);
35af33a6ddSJed Brown }
36af33a6ddSJed Brown 
379371c9d4SSatish Balay PetscErrorCode CharacteristicDestroy(Characteristic *c) {
38af33a6ddSJed Brown   PetscFunctionBegin;
396bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
406bf464f9SBarry Smith   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
416bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0);
42af33a6ddSJed Brown 
43*48a46eb9SPierre Jolivet   if ((*c)->ops->destroy) PetscCall((*(*c)->ops->destroy)((*c)));
449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&(*c)->itemType));
459566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queue));
469566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queueLocal));
479566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queueRemote));
489566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->neighbors));
499566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->needCount));
509566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->localOffsets));
519566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->fillCount));
529566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->remoteOffsets));
539566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->request));
549566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->status));
559566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(c));
56af33a6ddSJed Brown   PetscFunctionReturn(0);
57af33a6ddSJed Brown }
58af33a6ddSJed Brown 
599371c9d4SSatish Balay PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c) {
60af33a6ddSJed Brown   Characteristic newC;
61af33a6ddSJed Brown 
62af33a6ddSJed Brown   PetscFunctionBegin;
63af33a6ddSJed Brown   PetscValidPointer(c, 2);
640298fd71SBarry Smith   *c = NULL;
659566063dSJacob Faibussowitsch   PetscCall(CharacteristicInitializePackage());
66af33a6ddSJed Brown 
679566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "Characteristic", comm, CharacteristicDestroy, CharacteristicView));
68af33a6ddSJed Brown   *c = newC;
69af33a6ddSJed Brown 
70af33a6ddSJed Brown   newC->structured          = PETSC_TRUE;
71af33a6ddSJed Brown   newC->numIds              = 0;
720298fd71SBarry Smith   newC->velocityDA          = NULL;
730298fd71SBarry Smith   newC->velocity            = NULL;
740298fd71SBarry Smith   newC->velocityOld         = NULL;
75af33a6ddSJed Brown   newC->numVelocityComp     = 0;
760298fd71SBarry Smith   newC->velocityComp        = NULL;
770298fd71SBarry Smith   newC->velocityInterp      = NULL;
780298fd71SBarry Smith   newC->velocityInterpLocal = NULL;
790298fd71SBarry Smith   newC->velocityCtx         = NULL;
800298fd71SBarry Smith   newC->fieldDA             = NULL;
810298fd71SBarry Smith   newC->field               = NULL;
82af33a6ddSJed Brown   newC->numFieldComp        = 0;
830298fd71SBarry Smith   newC->fieldComp           = NULL;
840298fd71SBarry Smith   newC->fieldInterp         = NULL;
850298fd71SBarry Smith   newC->fieldInterpLocal    = NULL;
860298fd71SBarry Smith   newC->fieldCtx            = NULL;
870298fd71SBarry Smith   newC->itemType            = 0;
880298fd71SBarry Smith   newC->queue               = NULL;
89af33a6ddSJed Brown   newC->queueSize           = 0;
90af33a6ddSJed Brown   newC->queueMax            = 0;
910298fd71SBarry Smith   newC->queueLocal          = NULL;
92af33a6ddSJed Brown   newC->queueLocalSize      = 0;
93af33a6ddSJed Brown   newC->queueLocalMax       = 0;
940298fd71SBarry Smith   newC->queueRemote         = NULL;
95af33a6ddSJed Brown   newC->queueRemoteSize     = 0;
96af33a6ddSJed Brown   newC->queueRemoteMax      = 0;
97af33a6ddSJed Brown   newC->numNeighbors        = 0;
980298fd71SBarry Smith   newC->neighbors           = NULL;
990298fd71SBarry Smith   newC->needCount           = NULL;
1000298fd71SBarry Smith   newC->localOffsets        = NULL;
1010298fd71SBarry Smith   newC->fillCount           = NULL;
1020298fd71SBarry Smith   newC->remoteOffsets       = NULL;
1030298fd71SBarry Smith   newC->request             = NULL;
1040298fd71SBarry Smith   newC->status              = NULL;
105af33a6ddSJed Brown   PetscFunctionReturn(0);
106af33a6ddSJed Brown }
107af33a6ddSJed Brown 
108af33a6ddSJed Brown /*@C
109af33a6ddSJed Brown    CharacteristicSetType - Builds Characteristic for a particular solver.
110af33a6ddSJed Brown 
111af33a6ddSJed Brown    Logically Collective on Characteristic
112af33a6ddSJed Brown 
113af33a6ddSJed Brown    Input Parameters:
114af33a6ddSJed Brown +  c    - the method of characteristics context
115af33a6ddSJed Brown -  type - a known method
116af33a6ddSJed Brown 
117af33a6ddSJed Brown    Options Database Key:
118af33a6ddSJed Brown .  -characteristic_type <method> - Sets the method; use -help for a list
119af33a6ddSJed Brown     of available methods
120af33a6ddSJed Brown 
121af33a6ddSJed Brown    Notes:
12230a76a96SBarry Smith    See "include/petsccharacteristic.h" for available methods
123af33a6ddSJed Brown 
124af33a6ddSJed Brown   Normally, it is best to use the CharacteristicSetFromOptions() command and
125af33a6ddSJed Brown   then set the Characteristic type from the options database rather than by using
126af33a6ddSJed Brown   this routine.  Using the options database provides the user with
127af33a6ddSJed Brown   maximum flexibility in evaluating the many different Krylov methods.
128af33a6ddSJed Brown   The CharacteristicSetType() routine is provided for those situations where it
129af33a6ddSJed Brown   is necessary to set the iterative solver independently of the command
130af33a6ddSJed Brown   line or options database.  This might be the case, for example, when
131af33a6ddSJed Brown   the choice of iterative solver changes during the execution of the
132af33a6ddSJed Brown   program, and the user's application is taking responsibility for
133af33a6ddSJed Brown   choosing the appropriate method.  In other words, this routine is
134af33a6ddSJed Brown   not for beginners.
135af33a6ddSJed Brown 
136af33a6ddSJed Brown   Level: intermediate
137af33a6ddSJed Brown 
138db781477SPatrick Sanan .seealso: `CharacteristicType`
139af33a6ddSJed Brown 
140af33a6ddSJed Brown @*/
1419371c9d4SSatish Balay PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type) {
142af33a6ddSJed Brown   PetscBool match;
1435f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(Characteristic);
144af33a6ddSJed Brown 
145af33a6ddSJed Brown   PetscFunctionBegin;
146af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
147af33a6ddSJed Brown   PetscValidCharPointer(type, 2);
148af33a6ddSJed Brown 
1499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)c, type, &match));
150af33a6ddSJed Brown   if (match) PetscFunctionReturn(0);
151af33a6ddSJed Brown 
152af33a6ddSJed Brown   if (c->data) {
153af33a6ddSJed Brown     /* destroy the old private Characteristic context */
154dbbe0bcdSBarry Smith     PetscUseTypeMethod(c, destroy);
1550298fd71SBarry Smith     c->ops->destroy = NULL;
1562120983fSLisandro Dalcin     c->data         = NULL;
157af33a6ddSJed Brown   }
158af33a6ddSJed Brown 
1599566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(CharacteristicList, type, &r));
1603c633725SBarry Smith   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
161af33a6ddSJed Brown   c->setupcalled = 0;
1629566063dSJacob Faibussowitsch   PetscCall((*r)(c));
1639566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)c, type));
164af33a6ddSJed Brown   PetscFunctionReturn(0);
165af33a6ddSJed Brown }
166af33a6ddSJed Brown 
167af33a6ddSJed Brown /*@
168af33a6ddSJed Brown    CharacteristicSetUp - Sets up the internal data structures for the
169af33a6ddSJed Brown    later use of an iterative solver.
170af33a6ddSJed Brown 
171af33a6ddSJed Brown    Collective on Characteristic
172af33a6ddSJed Brown 
173af33a6ddSJed Brown    Input Parameter:
174af33a6ddSJed Brown .  ksp   - iterative context obtained from CharacteristicCreate()
175af33a6ddSJed Brown 
176af33a6ddSJed Brown    Level: developer
177af33a6ddSJed Brown 
178db781477SPatrick Sanan .seealso: `CharacteristicCreate()`, `CharacteristicSolve()`, `CharacteristicDestroy()`
179af33a6ddSJed Brown @*/
1809371c9d4SSatish Balay PetscErrorCode CharacteristicSetUp(Characteristic c) {
181af33a6ddSJed Brown   PetscFunctionBegin;
182af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
183af33a6ddSJed Brown 
184*48a46eb9SPierre Jolivet   if (!((PetscObject)c)->type_name) PetscCall(CharacteristicSetType(c, CHARACTERISTICDA));
185af33a6ddSJed Brown 
186af33a6ddSJed Brown   if (c->setupcalled == 2) PetscFunctionReturn(0);
187af33a6ddSJed Brown 
1889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
1899371c9d4SSatish Balay   if (!c->setupcalled) { PetscUseTypeMethod(c, setup); }
1909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
191af33a6ddSJed Brown   c->setupcalled = 2;
192af33a6ddSJed Brown   PetscFunctionReturn(0);
193af33a6ddSJed Brown }
194af33a6ddSJed Brown 
195af33a6ddSJed Brown /*@C
1961c84c290SBarry Smith   CharacteristicRegister -  Adds a solver to the method of characteristics package.
1971c84c290SBarry Smith 
1981c84c290SBarry Smith    Not Collective
1991c84c290SBarry Smith 
2001c84c290SBarry Smith    Input Parameters:
2011c84c290SBarry Smith +  name_solver - name of a new user-defined solver
2021c84c290SBarry Smith -  routine_create - routine to create method context
2031c84c290SBarry Smith 
2041c84c290SBarry Smith   Sample usage:
2051c84c290SBarry Smith .vb
206bdf89e91SBarry Smith     CharacteristicRegister("my_char", MyCharCreate);
2071c84c290SBarry Smith .ve
2081c84c290SBarry Smith 
2091c84c290SBarry Smith   Then, your Characteristic type can be chosen with the procedural interface via
2101c84c290SBarry Smith .vb
2111c84c290SBarry Smith     CharacteristicCreate(MPI_Comm, Characteristic* &char);
2121c84c290SBarry Smith     CharacteristicSetType(char,"my_char");
2131c84c290SBarry Smith .ve
2141c84c290SBarry Smith    or at runtime via the option
2151c84c290SBarry Smith .vb
2161c84c290SBarry Smith     -characteristic_type my_char
2171c84c290SBarry Smith .ve
2181c84c290SBarry Smith 
2191c84c290SBarry Smith    Notes:
2201c84c290SBarry Smith    CharacteristicRegister() may be called multiple times to add several user-defined solvers.
2211c84c290SBarry Smith 
222db781477SPatrick Sanan .seealso: `CharacteristicRegisterAll()`, `CharacteristicRegisterDestroy()`
223af33a6ddSJed Brown 
224af33a6ddSJed Brown   Level: advanced
225af33a6ddSJed Brown @*/
2269371c9d4SSatish Balay PetscErrorCode CharacteristicRegister(const char sname[], PetscErrorCode (*function)(Characteristic)) {
227af33a6ddSJed Brown   PetscFunctionBegin;
2289566063dSJacob Faibussowitsch   PetscCall(CharacteristicInitializePackage());
2299566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&CharacteristicList, sname, function));
230af33a6ddSJed Brown   PetscFunctionReturn(0);
231af33a6ddSJed Brown }
232af33a6ddSJed Brown 
2339371c9d4SSatish Balay PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx) {
234af33a6ddSJed Brown   PetscFunctionBegin;
235af33a6ddSJed Brown   c->velocityDA      = da;
236af33a6ddSJed Brown   c->velocity        = v;
237af33a6ddSJed Brown   c->velocityOld     = vOld;
238af33a6ddSJed Brown   c->numVelocityComp = numComponents;
239af33a6ddSJed Brown   c->velocityComp    = components;
240af33a6ddSJed Brown   c->velocityInterp  = interp;
241af33a6ddSJed Brown   c->velocityCtx     = ctx;
242af33a6ddSJed Brown   PetscFunctionReturn(0);
243af33a6ddSJed Brown }
244af33a6ddSJed Brown 
2459371c9d4SSatish Balay PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx) {
246af33a6ddSJed Brown   PetscFunctionBegin;
247af33a6ddSJed Brown   c->velocityDA          = da;
248af33a6ddSJed Brown   c->velocity            = v;
249af33a6ddSJed Brown   c->velocityOld         = vOld;
250af33a6ddSJed Brown   c->numVelocityComp     = numComponents;
251af33a6ddSJed Brown   c->velocityComp        = components;
252af33a6ddSJed Brown   c->velocityInterpLocal = interp;
253af33a6ddSJed Brown   c->velocityCtx         = ctx;
254af33a6ddSJed Brown   PetscFunctionReturn(0);
255af33a6ddSJed Brown }
256af33a6ddSJed Brown 
2579371c9d4SSatish Balay PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx) {
258af33a6ddSJed Brown   PetscFunctionBegin;
259af33a6ddSJed Brown #if 0
2603c633725SBarry 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.");
261af33a6ddSJed Brown #endif
262af33a6ddSJed Brown   c->fieldDA      = da;
263af33a6ddSJed Brown   c->field        = v;
264af33a6ddSJed Brown   c->numFieldComp = numComponents;
265af33a6ddSJed Brown   c->fieldComp    = components;
266af33a6ddSJed Brown   c->fieldInterp  = interp;
267af33a6ddSJed Brown   c->fieldCtx     = ctx;
268af33a6ddSJed Brown   PetscFunctionReturn(0);
269af33a6ddSJed Brown }
270af33a6ddSJed Brown 
2719371c9d4SSatish Balay PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx) {
272af33a6ddSJed Brown   PetscFunctionBegin;
273af33a6ddSJed Brown #if 0
2743c633725SBarry 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.");
275af33a6ddSJed Brown #endif
276af33a6ddSJed Brown   c->fieldDA          = da;
277af33a6ddSJed Brown   c->field            = v;
278af33a6ddSJed Brown   c->numFieldComp     = numComponents;
279af33a6ddSJed Brown   c->fieldComp        = components;
280af33a6ddSJed Brown   c->fieldInterpLocal = interp;
281af33a6ddSJed Brown   c->fieldCtx         = ctx;
282af33a6ddSJed Brown   PetscFunctionReturn(0);
283af33a6ddSJed Brown }
284af33a6ddSJed Brown 
2859371c9d4SSatish Balay PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution) {
286af33a6ddSJed Brown   CharacteristicPointDA2D Qi;
287af33a6ddSJed Brown   DM                      da = c->velocityDA;
288af33a6ddSJed Brown   Vec                     velocityLocal, velocityLocalOld;
289af33a6ddSJed Brown   Vec                     fieldLocal;
290af33a6ddSJed Brown   DMDALocalInfo           info;
291af33a6ddSJed Brown   PetscScalar           **solArray;
292af33a6ddSJed Brown   void                   *velocityArray;
293af33a6ddSJed Brown   void                   *velocityArrayOld;
294af33a6ddSJed Brown   void                   *fieldArray;
2951cc8b206SBarry Smith   PetscScalar            *interpIndices;
2961cc8b206SBarry Smith   PetscScalar            *velocityValues, *velocityValuesOld;
2971cc8b206SBarry Smith   PetscScalar            *fieldValues;
298af33a6ddSJed Brown   PetscMPIInt             rank;
299af33a6ddSJed Brown   PetscInt                dim;
300af33a6ddSJed Brown   PetscMPIInt             neighbors[9];
301af33a6ddSJed Brown   PetscInt                dof;
302af33a6ddSJed Brown   PetscInt                gx, gy;
303af33a6ddSJed Brown   PetscInt                n, is, ie, js, je, comp;
304af33a6ddSJed Brown 
305af33a6ddSJed Brown   PetscFunctionBegin;
306af33a6ddSJed Brown   c->queueSize = 0;
3079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
3089566063dSJacob Faibussowitsch   PetscCall(DMDAGetNeighborsRank(da, neighbors));
3099566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetNeighbors(c, 9, neighbors));
3109566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetUp(c));
311af33a6ddSJed Brown   /* global and local grid info */
3129566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, &gx, &gy, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3139566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
3149371c9d4SSatish Balay   is = info.xs;
3159371c9d4SSatish Balay   ie = info.xs + info.xm;
3169371c9d4SSatish Balay   js = info.ys;
3179371c9d4SSatish Balay   je = info.ys + info.ym;
318af33a6ddSJed Brown   /* Allocation */
3199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim, &interpIndices));
3209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValues));
3219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValuesOld));
3229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numFieldComp, &fieldValues));
3239566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
324af33a6ddSJed Brown 
325af33a6ddSJed Brown   /* -----------------------------------------------------------------------
326af33a6ddSJed Brown      PART 1, AT t-dt/2
327af33a6ddSJed Brown      -----------------------------------------------------------------------*/
3289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
329af33a6ddSJed Brown   /* GET POSITION AT HALF TIME IN THE PAST */
330af33a6ddSJed Brown   if (c->velocityInterpLocal) {
3319566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocal));
3329566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocalOld));
3339566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3349566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3359566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3369566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3379566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray));
3389566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
339af33a6ddSJed Brown   }
3409566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n"));
341af33a6ddSJed Brown   for (Qi.j = js; Qi.j < je; Qi.j++) {
342af33a6ddSJed Brown     for (Qi.i = is; Qi.i < ie; Qi.i++) {
343af33a6ddSJed Brown       interpIndices[0] = Qi.i;
344af33a6ddSJed Brown       interpIndices[1] = Qi.j;
3459566063dSJacob Faibussowitsch       if (c->velocityInterpLocal) PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3469566063dSJacob Faibussowitsch       else PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
347af33a6ddSJed Brown       Qi.x = Qi.i - velocityValues[0] * dt / 2.0;
348af33a6ddSJed Brown       Qi.y = Qi.j - velocityValues[1] * dt / 2.0;
349af33a6ddSJed Brown 
350af33a6ddSJed Brown       /* Determine whether the position at t - dt/2 is local */
351af33a6ddSJed Brown       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
352af33a6ddSJed Brown 
353af33a6ddSJed Brown       /* Check for Periodic boundaries and move all periodic points back onto the domain */
3549566063dSJacob Faibussowitsch       PetscCall(DMDAMapCoordsToPeriodicDomain(da, &(Qi.x), &(Qi.y)));
3559566063dSJacob Faibussowitsch       PetscCall(CharacteristicAddPoint(c, &Qi));
356af33a6ddSJed Brown     }
357af33a6ddSJed Brown   }
3589566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
359af33a6ddSJed Brown 
3609566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
3619566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
3629566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
363af33a6ddSJed Brown 
3649566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
365af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (local values) */
3669566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n"));
367af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
368af33a6ddSJed Brown     Qi = c->queue[n];
369af33a6ddSJed Brown     if (c->neighbors[Qi.proc] == rank) {
370af33a6ddSJed Brown       interpIndices[0] = Qi.x;
371af33a6ddSJed Brown       interpIndices[1] = Qi.y;
372af33a6ddSJed Brown       if (c->velocityInterpLocal) {
3739566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3749566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
375af33a6ddSJed Brown       } else {
3769566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3779566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
378af33a6ddSJed Brown       }
379af33a6ddSJed Brown       Qi.x = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
380af33a6ddSJed Brown       Qi.y = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
381af33a6ddSJed Brown     }
382af33a6ddSJed Brown     c->queue[n] = Qi;
383af33a6ddSJed Brown   }
3849566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
385af33a6ddSJed Brown 
3869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
3879566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
3889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
389af33a6ddSJed Brown 
390af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (fill remote requests) */
3919566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
39263a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote velocities at t_{n - 1/2}\n", c->queueRemoteSize));
393af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
394af33a6ddSJed Brown     Qi               = c->queueRemote[n];
395af33a6ddSJed Brown     interpIndices[0] = Qi.x;
396af33a6ddSJed Brown     interpIndices[1] = Qi.y;
397af33a6ddSJed Brown     if (c->velocityInterpLocal) {
3989566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3999566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
400af33a6ddSJed Brown     } else {
4019566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4029566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
403af33a6ddSJed Brown     }
404af33a6ddSJed Brown     Qi.x              = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
405af33a6ddSJed Brown     Qi.y              = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
406af33a6ddSJed Brown     c->queueRemote[n] = Qi;
407af33a6ddSJed Brown   }
4089566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
4099566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
4109566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
4119566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
412af33a6ddSJed Brown   if (c->velocityInterpLocal) {
4139566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray));
4149566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
4159566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocal));
4169566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocalOld));
417af33a6ddSJed Brown   }
4189566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
419af33a6ddSJed Brown 
420af33a6ddSJed Brown   /* -----------------------------------------------------------------------
421af33a6ddSJed Brown      PART 2, AT t-dt
422af33a6ddSJed Brown      -----------------------------------------------------------------------*/
423af33a6ddSJed Brown 
424af33a6ddSJed Brown   /* GET POSITION AT t_n (local values) */
4259566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
4269566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n}\n"));
427af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
428af33a6ddSJed Brown     Qi   = c->queue[n];
429af33a6ddSJed Brown     Qi.x = Qi.i - Qi.x * dt;
430af33a6ddSJed Brown     Qi.y = Qi.j - Qi.y * dt;
431af33a6ddSJed Brown 
432af33a6ddSJed Brown     /* Determine whether the position at t-dt is local */
433af33a6ddSJed Brown     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
434af33a6ddSJed Brown 
435af33a6ddSJed Brown     /* Check for Periodic boundaries and move all periodic points back onto the domain */
4369566063dSJacob Faibussowitsch     PetscCall(DMDAMapCoordsToPeriodicDomain(da, &(Qi.x), &(Qi.y)));
437af33a6ddSJed Brown 
438af33a6ddSJed Brown     c->queue[n] = Qi;
439af33a6ddSJed Brown   }
4409566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
441af33a6ddSJed Brown 
4429566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4439566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
4449566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
445af33a6ddSJed Brown 
446af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
4479566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
448af33a6ddSJed Brown   if (c->fieldInterpLocal) {
4499566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->fieldDA, &fieldLocal));
4509566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4519566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4529566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray));
453af33a6ddSJed Brown   }
4549566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local field at t_{n}\n"));
455af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
456af33a6ddSJed Brown     if (c->neighbors[c->queue[n].proc] == rank) {
457af33a6ddSJed Brown       interpIndices[0] = c->queue[n].x;
458af33a6ddSJed Brown       interpIndices[1] = c->queue[n].y;
4599566063dSJacob Faibussowitsch       if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4609566063dSJacob Faibussowitsch       else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
461bbd56ea5SKarl Rupp       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
462af33a6ddSJed Brown     }
463af33a6ddSJed Brown   }
4649566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
465af33a6ddSJed Brown 
4669566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4679566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
4689566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
469af33a6ddSJed Brown 
470af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
4719566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
47263a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote field points at t_{n}\n", c->queueRemoteSize));
473af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
474af33a6ddSJed Brown     interpIndices[0] = c->queueRemote[n].x;
475af33a6ddSJed Brown     interpIndices[1] = c->queueRemote[n].y;
476af33a6ddSJed Brown 
477af33a6ddSJed Brown     /* for debugging purposes */
478af33a6ddSJed Brown     if (1) { /* hacked bounds test...let's do better */
4799371c9d4SSatish Balay       PetscScalar im = interpIndices[0];
4809371c9d4SSatish Balay       PetscScalar jm = interpIndices[1];
481af33a6ddSJed Brown 
48263a3b9bcSJacob 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));
483af33a6ddSJed Brown     }
484af33a6ddSJed Brown 
4859566063dSJacob Faibussowitsch     if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4869566063dSJacob Faibussowitsch     else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
487bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
488af33a6ddSJed Brown   }
4899566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
490af33a6ddSJed Brown 
4919566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4929566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
4939566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
494af33a6ddSJed Brown   if (c->fieldInterpLocal) {
4959566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray));
4969566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->fieldDA, &fieldLocal));
497af33a6ddSJed Brown   }
4989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
499af33a6ddSJed Brown 
500af33a6ddSJed Brown   /* Return field of characteristics at t_n-1 */
5019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5029566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(c->fieldDA, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5039566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(c->fieldDA, solution, &solArray));
504af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
505af33a6ddSJed Brown     Qi = c->queue[n];
506bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i * dof + c->fieldComp[comp]] = Qi.field[comp];
507af33a6ddSJed Brown   }
5089566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(c->fieldDA, solution, &solArray));
5099566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
511af33a6ddSJed Brown 
512af33a6ddSJed Brown   /* Cleanup */
5139566063dSJacob Faibussowitsch   PetscCall(PetscFree(interpIndices));
5149566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValues));
5159566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValuesOld));
5169566063dSJacob Faibussowitsch   PetscCall(PetscFree(fieldValues));
517af33a6ddSJed Brown   PetscFunctionReturn(0);
518af33a6ddSJed Brown }
519af33a6ddSJed Brown 
5209371c9d4SSatish Balay PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[]) {
521af33a6ddSJed Brown   PetscFunctionBegin;
522af33a6ddSJed Brown   c->numNeighbors = numNeighbors;
5239566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->neighbors));
5249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numNeighbors, &c->neighbors));
5259566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(c->neighbors, neighbors, numNeighbors));
526af33a6ddSJed Brown   PetscFunctionReturn(0);
527af33a6ddSJed Brown }
528af33a6ddSJed Brown 
5299371c9d4SSatish Balay PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point) {
530af33a6ddSJed Brown   PetscFunctionBegin;
53163a3b9bcSJacob Faibussowitsch   PetscCheck(c->queueSize < c->queueMax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exceeded maximum queue size %" PetscInt_FMT, c->queueMax);
532af33a6ddSJed Brown   c->queue[c->queueSize++] = *point;
533af33a6ddSJed Brown   PetscFunctionReturn(0);
534af33a6ddSJed Brown }
535af33a6ddSJed Brown 
5369371c9d4SSatish Balay int CharacteristicSendCoordinatesBegin(Characteristic c) {
537af33a6ddSJed Brown   PetscMPIInt rank, tag = 121;
538af33a6ddSJed Brown   PetscInt    i, n;
539af33a6ddSJed Brown 
540af33a6ddSJed Brown   PetscFunctionBegin;
5419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
5429566063dSJacob Faibussowitsch   PetscCall(CharacteristicHeapSort(c, c->queue, c->queueSize));
5439566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(c->needCount, c->numNeighbors));
544bbd56ea5SKarl Rupp   for (i = 0; i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
545af33a6ddSJed Brown   c->fillCount[0] = 0;
546*48a46eb9SPierre 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])));
547*48a46eb9SPierre Jolivet   for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
5489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
549af33a6ddSJed Brown   /* Initialize the remote queue */
550af33a6ddSJed Brown   c->queueLocalMax = c->localOffsets[0] = 0;
551af33a6ddSJed Brown   c->queueRemoteMax = c->remoteOffsets[0] = 0;
552af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
553af33a6ddSJed Brown     c->remoteOffsets[n] = c->queueRemoteMax;
554af33a6ddSJed Brown     c->queueRemoteMax += c->fillCount[n];
555af33a6ddSJed Brown     c->localOffsets[n] = c->queueLocalMax;
556af33a6ddSJed Brown     c->queueLocalMax += c->needCount[n];
557af33a6ddSJed Brown   }
558af33a6ddSJed Brown   /* HACK BEGIN */
559bbd56ea5SKarl Rupp   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
560af33a6ddSJed Brown   c->needCount[0] = 0;
561af33a6ddSJed Brown   /* HACK END */
562af33a6ddSJed Brown   if (c->queueRemoteMax) {
5639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(c->queueRemoteMax, &c->queueRemote));
5640298fd71SBarry Smith   } else c->queueRemote = NULL;
565af33a6ddSJed Brown   c->queueRemoteSize = c->queueRemoteMax;
566af33a6ddSJed Brown 
567af33a6ddSJed Brown   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
568af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
56963a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Receiving %" PetscInt_FMT " requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]));
5709566063dSJacob 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])));
571af33a6ddSJed Brown   }
572af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
57363a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Sending %" PetscInt_FMT " requests for values from proc %d\n", c->needCount[n], c->neighbors[n]));
5749566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
575af33a6ddSJed Brown   }
576af33a6ddSJed Brown   PetscFunctionReturn(0);
577af33a6ddSJed Brown }
578af33a6ddSJed Brown 
5799371c9d4SSatish Balay PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c) {
580af33a6ddSJed Brown #if 0
581af33a6ddSJed Brown   PetscMPIInt rank;
582af33a6ddSJed Brown   PetscInt    n;
583af33a6ddSJed Brown #endif
584af33a6ddSJed Brown 
585af33a6ddSJed Brown   PetscFunctionBegin;
5869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
587af33a6ddSJed Brown #if 0
5889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
589af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
5903c633725SBarry 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);
591af33a6ddSJed Brown   }
592af33a6ddSJed Brown #endif
593af33a6ddSJed Brown   PetscFunctionReturn(0);
594af33a6ddSJed Brown }
595af33a6ddSJed Brown 
5969371c9d4SSatish Balay PetscErrorCode CharacteristicGetValuesBegin(Characteristic c) {
597af33a6ddSJed Brown   PetscMPIInt tag = 121;
598af33a6ddSJed Brown   PetscInt    n;
599af33a6ddSJed Brown 
600af33a6ddSJed Brown   PetscFunctionBegin;
601af33a6ddSJed Brown   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
602*48a46eb9SPierre 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])));
603*48a46eb9SPierre 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)));
604af33a6ddSJed Brown   PetscFunctionReturn(0);
605af33a6ddSJed Brown }
606af33a6ddSJed Brown 
6079371c9d4SSatish Balay PetscErrorCode CharacteristicGetValuesEnd(Characteristic c) {
608af33a6ddSJed Brown   PetscFunctionBegin;
6099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
610af33a6ddSJed Brown   /* Free queue of requests from other procs */
6119566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->queueRemote));
612af33a6ddSJed Brown   PetscFunctionReturn(0);
613af33a6ddSJed Brown }
614af33a6ddSJed Brown 
615af33a6ddSJed Brown /*---------------------------------------------------------------------*/
616af33a6ddSJed Brown /*
617af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
618af33a6ddSJed Brown */
6190b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
620af33a6ddSJed Brown /*---------------------------------------------------------------------*/
621af33a6ddSJed Brown {
622af33a6ddSJed Brown   CharacteristicPointDA2D temp;
6230b8f38c8SBarry Smith   PetscInt                n;
624af33a6ddSJed Brown 
6250b8f38c8SBarry Smith   PetscFunctionBegin;
626af33a6ddSJed Brown   if (0) { /* Check the order of the queue before sorting */
6279566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Before Heap sort\n"));
628*48a46eb9SPierre Jolivet     for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
629af33a6ddSJed Brown   }
630af33a6ddSJed Brown 
631af33a6ddSJed Brown   /* SORTING PHASE */
6329371c9d4SSatish Balay   for (n = (size / 2) - 1; n >= 0; n--) { PetscCall(CharacteristicSiftDown(c, queue, n, size - 1)); /* Rich had size-1 here, Matt had size*/ }
633af33a6ddSJed Brown   for (n = size - 1; n >= 1; n--) {
634af33a6ddSJed Brown     temp     = queue[0];
635af33a6ddSJed Brown     queue[0] = queue[n];
636af33a6ddSJed Brown     queue[n] = temp;
6379566063dSJacob Faibussowitsch     PetscCall(CharacteristicSiftDown(c, queue, 0, n - 1));
638af33a6ddSJed Brown   }
639af33a6ddSJed Brown   if (0) { /* Check the order of the queue after sorting */
6409566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Avter  Heap sort\n"));
641*48a46eb9SPierre Jolivet     for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
642af33a6ddSJed Brown   }
6430b8f38c8SBarry Smith   PetscFunctionReturn(0);
644af33a6ddSJed Brown }
645af33a6ddSJed Brown 
646af33a6ddSJed Brown /*---------------------------------------------------------------------*/
647af33a6ddSJed Brown /*
648af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
649af33a6ddSJed Brown */
6500b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
651af33a6ddSJed Brown /*---------------------------------------------------------------------*/
652af33a6ddSJed Brown {
653af33a6ddSJed Brown   PetscBool               done = PETSC_FALSE;
654af33a6ddSJed Brown   PetscInt                maxChild;
655af33a6ddSJed Brown   CharacteristicPointDA2D temp;
656af33a6ddSJed Brown 
657af33a6ddSJed Brown   PetscFunctionBegin;
658af33a6ddSJed Brown   while ((root * 2 <= bottom) && (!done)) {
659af33a6ddSJed Brown     if (root * 2 == bottom) maxChild = root * 2;
660af33a6ddSJed Brown     else if (queue[root * 2].proc > queue[root * 2 + 1].proc) maxChild = root * 2;
661af33a6ddSJed Brown     else maxChild = root * 2 + 1;
662af33a6ddSJed Brown 
663af33a6ddSJed Brown     if (queue[root].proc < queue[maxChild].proc) {
664af33a6ddSJed Brown       temp            = queue[root];
665af33a6ddSJed Brown       queue[root]     = queue[maxChild];
666af33a6ddSJed Brown       queue[maxChild] = temp;
667af33a6ddSJed Brown       root            = maxChild;
668db4deed7SKarl Rupp     } else done = PETSC_TRUE;
669af33a6ddSJed Brown   }
670af33a6ddSJed Brown   PetscFunctionReturn(0);
671af33a6ddSJed Brown }
672af33a6ddSJed Brown 
673af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
6749371c9d4SSatish Balay PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[]) {
675bff4a2f0SMatthew G. Knepley   DMBoundaryType bx, by;
676af33a6ddSJed Brown   PetscBool      IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
677af33a6ddSJed Brown   MPI_Comm       comm;
678af33a6ddSJed Brown   PetscMPIInt    rank;
679af33a6ddSJed Brown   PetscInt     **procs, pi, pj, pim, pip, pjm, pjp, PI, PJ;
680af33a6ddSJed Brown 
681af33a6ddSJed Brown   PetscFunctionBegin;
6829566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
6839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
6849566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, &PI, &PJ, NULL, NULL, NULL, &bx, &by, NULL, NULL));
685af33a6ddSJed Brown 
686bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
687bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
688af33a6ddSJed Brown 
689af33a6ddSJed Brown   neighbors[0] = rank;
690af33a6ddSJed Brown   rank         = 0;
6919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(PJ, &procs));
692af33a6ddSJed Brown   for (pj = 0; pj < PJ; pj++) {
6939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PI, &(procs[pj])));
694af33a6ddSJed Brown     for (pi = 0; pi < PI; pi++) {
695af33a6ddSJed Brown       procs[pj][pi] = rank;
696af33a6ddSJed Brown       rank++;
697af33a6ddSJed Brown     }
698af33a6ddSJed Brown   }
699af33a6ddSJed Brown 
700af33a6ddSJed Brown   pi  = neighbors[0] % PI;
701af33a6ddSJed Brown   pj  = neighbors[0] / PI;
7029371c9d4SSatish Balay   pim = pi - 1;
7039371c9d4SSatish Balay   if (pim < 0) pim = PI - 1;
704af33a6ddSJed Brown   pip = (pi + 1) % PI;
7059371c9d4SSatish Balay   pjm = pj - 1;
7069371c9d4SSatish Balay   if (pjm < 0) pjm = PJ - 1;
707af33a6ddSJed Brown   pjp = (pj + 1) % PJ;
708af33a6ddSJed Brown 
709af33a6ddSJed Brown   neighbors[1] = procs[pj][pim];
710af33a6ddSJed Brown   neighbors[2] = procs[pjp][pim];
711af33a6ddSJed Brown   neighbors[3] = procs[pjp][pi];
712af33a6ddSJed Brown   neighbors[4] = procs[pjp][pip];
713af33a6ddSJed Brown   neighbors[5] = procs[pj][pip];
714af33a6ddSJed Brown   neighbors[6] = procs[pjm][pip];
715af33a6ddSJed Brown   neighbors[7] = procs[pjm][pi];
716af33a6ddSJed Brown   neighbors[8] = procs[pjm][pim];
717af33a6ddSJed Brown 
718af33a6ddSJed Brown   if (!IPeriodic) {
719af33a6ddSJed Brown     if (pi == 0) neighbors[1] = neighbors[2] = neighbors[8] = neighbors[0];
720af33a6ddSJed Brown     if (pi == PI - 1) neighbors[4] = neighbors[5] = neighbors[6] = neighbors[0];
721af33a6ddSJed Brown   }
722af33a6ddSJed Brown 
723af33a6ddSJed Brown   if (!JPeriodic) {
724af33a6ddSJed Brown     if (pj == 0) neighbors[6] = neighbors[7] = neighbors[8] = neighbors[0];
725af33a6ddSJed Brown     if (pj == PJ - 1) neighbors[2] = neighbors[3] = neighbors[4] = neighbors[0];
726af33a6ddSJed Brown   }
727af33a6ddSJed Brown 
728*48a46eb9SPierre Jolivet   for (pj = 0; pj < PJ; pj++) PetscCall(PetscFree(procs[pj]));
7299566063dSJacob Faibussowitsch   PetscCall(PetscFree(procs));
730af33a6ddSJed Brown   PetscFunctionReturn(0);
731af33a6ddSJed Brown }
732af33a6ddSJed Brown 
733af33a6ddSJed Brown /*
734af33a6ddSJed Brown   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
735af33a6ddSJed Brown     2 | 3 | 4
736af33a6ddSJed Brown     __|___|__
737af33a6ddSJed Brown     1 | 0 | 5
738af33a6ddSJed Brown     __|___|__
739af33a6ddSJed Brown     8 | 7 | 6
740af33a6ddSJed Brown       |   |
741af33a6ddSJed Brown */
7429371c9d4SSatish Balay PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr) {
743af33a6ddSJed Brown   DMDALocalInfo info;
7441cc8b206SBarry Smith   PetscReal     is, ie, js, je;
745af33a6ddSJed Brown 
7469566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
7479371c9d4SSatish Balay   is = (PetscReal)info.xs - 0.5;
7489371c9d4SSatish Balay   ie = (PetscReal)info.xs + info.xm - 0.5;
7499371c9d4SSatish Balay   js = (PetscReal)info.ys - 0.5;
7509371c9d4SSatish Balay   je = (PetscReal)info.ys + info.ym - 0.5;
751af33a6ddSJed Brown 
752af33a6ddSJed Brown   if (ir >= is && ir <= ie) { /* center column */
753bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 0;
754bbd56ea5SKarl Rupp     else if (jr < js) return 7;
755bbd56ea5SKarl Rupp     else return 3;
756af33a6ddSJed Brown   } else if (ir < is) { /* left column */
757bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 1;
758bbd56ea5SKarl Rupp     else if (jr < js) return 8;
759bbd56ea5SKarl Rupp     else return 2;
760af33a6ddSJed Brown   } else { /* right column */
761bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 5;
762bbd56ea5SKarl Rupp     else if (jr < js) return 6;
763bbd56ea5SKarl Rupp     else return 4;
764af33a6ddSJed Brown   }
765af33a6ddSJed Brown }
766