xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1af33a6ddSJed Brown 
2af0996ceSBarry Smith #include <petsc/private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/
31e25c274SJed Brown #include <petscdmda.h>
4665c2dedSJed Brown #include <petscviewer.h>
5af33a6ddSJed Brown 
6af33a6ddSJed Brown PetscClassId  CHARACTERISTIC_CLASSID;
7af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate;
8af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange;
9af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange;
10af33a6ddSJed Brown /*
11af33a6ddSJed Brown    Contains the list of registered characteristic routines
12af33a6ddSJed Brown */
130298fd71SBarry Smith PetscFunctionList CharacteristicList              = NULL;
14140e18c1SBarry Smith PetscBool         CharacteristicRegisterAllCalled = PETSC_FALSE;
15af33a6ddSJed Brown 
16af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt []);
171cc8b206SBarry Smith PetscInt       DMDAGetNeighborRelative(DM, PetscReal, PetscReal);
18af33a6ddSJed Brown PetscErrorCode DMDAMapToPeriodicDomain(DM, PetscScalar []);
19af33a6ddSJed Brown 
200b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt);
210b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt);
22af33a6ddSJed Brown 
23af33a6ddSJed Brown PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
24af33a6ddSJed Brown {
25af33a6ddSJed Brown   PetscBool      iascii;
26af33a6ddSJed Brown 
27af33a6ddSJed Brown   PetscFunctionBegin;
28af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
29af33a6ddSJed Brown   if (!viewer) {
30*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer));
31af33a6ddSJed Brown   }
32af33a6ddSJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
33af33a6ddSJed Brown   PetscCheckSameComm(c, 1, viewer, 2);
34af33a6ddSJed Brown 
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
36af33a6ddSJed Brown   if (!iascii) {
37af33a6ddSJed Brown     if (c->ops->view) {
38*5f80ce2aSJacob Faibussowitsch       CHKERRQ((*c->ops->view)(c, viewer));
39af33a6ddSJed Brown     }
40af33a6ddSJed Brown   }
41af33a6ddSJed Brown   PetscFunctionReturn(0);
42af33a6ddSJed Brown }
43af33a6ddSJed Brown 
446bf464f9SBarry Smith PetscErrorCode CharacteristicDestroy(Characteristic *c)
45af33a6ddSJed Brown {
46af33a6ddSJed Brown   PetscFunctionBegin;
476bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
486bf464f9SBarry Smith   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
496bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0);
50af33a6ddSJed Brown 
516bf464f9SBarry Smith   if ((*c)->ops->destroy) {
52*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*(*c)->ops->destroy)((*c)));
53af33a6ddSJed Brown   }
54*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&(*c)->itemType));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->queue));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->queueLocal));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->queueRemote));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->neighbors));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->needCount));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->localOffsets));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->fillCount));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->remoteOffsets));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->request));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*c)->status));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHeaderDestroy(c));
66af33a6ddSJed Brown   PetscFunctionReturn(0);
67af33a6ddSJed Brown }
68af33a6ddSJed Brown 
69af33a6ddSJed Brown PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
70af33a6ddSJed Brown {
71af33a6ddSJed Brown   Characteristic newC;
72af33a6ddSJed Brown 
73af33a6ddSJed Brown   PetscFunctionBegin;
74af33a6ddSJed Brown   PetscValidPointer(c, 2);
750298fd71SBarry Smith   *c = NULL;
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicInitializePackage());
77af33a6ddSJed Brown 
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "Characteristic", comm, CharacteristicDestroy, CharacteristicView));
79af33a6ddSJed Brown   *c   = newC;
80af33a6ddSJed Brown 
81af33a6ddSJed Brown   newC->structured          = PETSC_TRUE;
82af33a6ddSJed Brown   newC->numIds              = 0;
830298fd71SBarry Smith   newC->velocityDA          = NULL;
840298fd71SBarry Smith   newC->velocity            = NULL;
850298fd71SBarry Smith   newC->velocityOld         = NULL;
86af33a6ddSJed Brown   newC->numVelocityComp     = 0;
870298fd71SBarry Smith   newC->velocityComp        = NULL;
880298fd71SBarry Smith   newC->velocityInterp      = NULL;
890298fd71SBarry Smith   newC->velocityInterpLocal = NULL;
900298fd71SBarry Smith   newC->velocityCtx         = NULL;
910298fd71SBarry Smith   newC->fieldDA             = NULL;
920298fd71SBarry Smith   newC->field               = NULL;
93af33a6ddSJed Brown   newC->numFieldComp        = 0;
940298fd71SBarry Smith   newC->fieldComp           = NULL;
950298fd71SBarry Smith   newC->fieldInterp         = NULL;
960298fd71SBarry Smith   newC->fieldInterpLocal    = NULL;
970298fd71SBarry Smith   newC->fieldCtx            = NULL;
980298fd71SBarry Smith   newC->itemType            = 0;
990298fd71SBarry Smith   newC->queue               = NULL;
100af33a6ddSJed Brown   newC->queueSize           = 0;
101af33a6ddSJed Brown   newC->queueMax            = 0;
1020298fd71SBarry Smith   newC->queueLocal          = NULL;
103af33a6ddSJed Brown   newC->queueLocalSize      = 0;
104af33a6ddSJed Brown   newC->queueLocalMax       = 0;
1050298fd71SBarry Smith   newC->queueRemote         = NULL;
106af33a6ddSJed Brown   newC->queueRemoteSize     = 0;
107af33a6ddSJed Brown   newC->queueRemoteMax      = 0;
108af33a6ddSJed Brown   newC->numNeighbors        = 0;
1090298fd71SBarry Smith   newC->neighbors           = NULL;
1100298fd71SBarry Smith   newC->needCount           = NULL;
1110298fd71SBarry Smith   newC->localOffsets        = NULL;
1120298fd71SBarry Smith   newC->fillCount           = NULL;
1130298fd71SBarry Smith   newC->remoteOffsets       = NULL;
1140298fd71SBarry Smith   newC->request             = NULL;
1150298fd71SBarry Smith   newC->status              = NULL;
116af33a6ddSJed Brown   PetscFunctionReturn(0);
117af33a6ddSJed Brown }
118af33a6ddSJed Brown 
119af33a6ddSJed Brown /*@C
120af33a6ddSJed Brown    CharacteristicSetType - Builds Characteristic for a particular solver.
121af33a6ddSJed Brown 
122af33a6ddSJed Brown    Logically Collective on Characteristic
123af33a6ddSJed Brown 
124af33a6ddSJed Brown    Input Parameters:
125af33a6ddSJed Brown +  c    - the method of characteristics context
126af33a6ddSJed Brown -  type - a known method
127af33a6ddSJed Brown 
128af33a6ddSJed Brown    Options Database Key:
129af33a6ddSJed Brown .  -characteristic_type <method> - Sets the method; use -help for a list
130af33a6ddSJed Brown     of available methods
131af33a6ddSJed Brown 
132af33a6ddSJed Brown    Notes:
13330a76a96SBarry Smith    See "include/petsccharacteristic.h" for available methods
134af33a6ddSJed Brown 
135af33a6ddSJed Brown   Normally, it is best to use the CharacteristicSetFromOptions() command and
136af33a6ddSJed Brown   then set the Characteristic type from the options database rather than by using
137af33a6ddSJed Brown   this routine.  Using the options database provides the user with
138af33a6ddSJed Brown   maximum flexibility in evaluating the many different Krylov methods.
139af33a6ddSJed Brown   The CharacteristicSetType() routine is provided for those situations where it
140af33a6ddSJed Brown   is necessary to set the iterative solver independently of the command
141af33a6ddSJed Brown   line or options database.  This might be the case, for example, when
142af33a6ddSJed Brown   the choice of iterative solver changes during the execution of the
143af33a6ddSJed Brown   program, and the user's application is taking responsibility for
144af33a6ddSJed Brown   choosing the appropriate method.  In other words, this routine is
145af33a6ddSJed Brown   not for beginners.
146af33a6ddSJed Brown 
147af33a6ddSJed Brown   Level: intermediate
148af33a6ddSJed Brown 
149af33a6ddSJed Brown .seealso: CharacteristicType
150af33a6ddSJed Brown 
151af33a6ddSJed Brown @*/
15219fd82e9SBarry Smith PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
153af33a6ddSJed Brown {
154af33a6ddSJed Brown   PetscBool      match;
155*5f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(Characteristic);
156af33a6ddSJed Brown 
157af33a6ddSJed Brown   PetscFunctionBegin;
158af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
159af33a6ddSJed Brown   PetscValidCharPointer(type, 2);
160af33a6ddSJed Brown 
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) c, type, &match));
162af33a6ddSJed Brown   if (match) PetscFunctionReturn(0);
163af33a6ddSJed Brown 
164af33a6ddSJed Brown   if (c->data) {
165af33a6ddSJed Brown     /* destroy the old private Characteristic context */
166*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*c->ops->destroy)(c));
1670298fd71SBarry Smith     c->ops->destroy = NULL;
1682120983fSLisandro Dalcin     c->data         = NULL;
169af33a6ddSJed Brown   }
170af33a6ddSJed Brown 
171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListFind(CharacteristicList,type,&r));
1723c633725SBarry Smith   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
173af33a6ddSJed Brown   c->setupcalled = 0;
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*r)(c));
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectChangeTypeName((PetscObject) c, type));
176af33a6ddSJed Brown   PetscFunctionReturn(0);
177af33a6ddSJed Brown }
178af33a6ddSJed Brown 
179af33a6ddSJed Brown /*@
180af33a6ddSJed Brown    CharacteristicSetUp - Sets up the internal data structures for the
181af33a6ddSJed Brown    later use of an iterative solver.
182af33a6ddSJed Brown 
183af33a6ddSJed Brown    Collective on Characteristic
184af33a6ddSJed Brown 
185af33a6ddSJed Brown    Input Parameter:
186af33a6ddSJed Brown .  ksp   - iterative context obtained from CharacteristicCreate()
187af33a6ddSJed Brown 
188af33a6ddSJed Brown    Level: developer
189af33a6ddSJed Brown 
190af33a6ddSJed Brown .seealso: CharacteristicCreate(), CharacteristicSolve(), CharacteristicDestroy()
191af33a6ddSJed Brown @*/
192af33a6ddSJed Brown PetscErrorCode CharacteristicSetUp(Characteristic c)
193af33a6ddSJed Brown {
194af33a6ddSJed Brown   PetscFunctionBegin;
195af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
196af33a6ddSJed Brown 
197af33a6ddSJed Brown   if (!((PetscObject)c)->type_name) {
198*5f80ce2aSJacob Faibussowitsch     CHKERRQ(CharacteristicSetType(c, CHARACTERISTICDA));
199af33a6ddSJed Brown   }
200af33a6ddSJed Brown 
201af33a6ddSJed Brown   if (c->setupcalled == 2) PetscFunctionReturn(0);
202af33a6ddSJed Brown 
203*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_SetUp,c,NULL,NULL,NULL));
204af33a6ddSJed Brown   if (!c->setupcalled) {
205*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*c->ops->setup)(c));
206af33a6ddSJed Brown   }
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_SetUp,c,NULL,NULL,NULL));
208af33a6ddSJed Brown   c->setupcalled = 2;
209af33a6ddSJed Brown   PetscFunctionReturn(0);
210af33a6ddSJed Brown }
211af33a6ddSJed Brown 
212af33a6ddSJed Brown /*@C
2131c84c290SBarry Smith   CharacteristicRegister -  Adds a solver to the method of characteristics package.
2141c84c290SBarry Smith 
2151c84c290SBarry Smith    Not Collective
2161c84c290SBarry Smith 
2171c84c290SBarry Smith    Input Parameters:
2181c84c290SBarry Smith +  name_solver - name of a new user-defined solver
2191c84c290SBarry Smith -  routine_create - routine to create method context
2201c84c290SBarry Smith 
2211c84c290SBarry Smith   Sample usage:
2221c84c290SBarry Smith .vb
223bdf89e91SBarry Smith     CharacteristicRegister("my_char", MyCharCreate);
2241c84c290SBarry Smith .ve
2251c84c290SBarry Smith 
2261c84c290SBarry Smith   Then, your Characteristic type can be chosen with the procedural interface via
2271c84c290SBarry Smith .vb
2281c84c290SBarry Smith     CharacteristicCreate(MPI_Comm, Characteristic* &char);
2291c84c290SBarry Smith     CharacteristicSetType(char,"my_char");
2301c84c290SBarry Smith .ve
2311c84c290SBarry Smith    or at runtime via the option
2321c84c290SBarry Smith .vb
2331c84c290SBarry Smith     -characteristic_type my_char
2341c84c290SBarry Smith .ve
2351c84c290SBarry Smith 
2361c84c290SBarry Smith    Notes:
2371c84c290SBarry Smith    CharacteristicRegister() may be called multiple times to add several user-defined solvers.
2381c84c290SBarry Smith 
2391c84c290SBarry Smith .seealso: CharacteristicRegisterAll(), CharacteristicRegisterDestroy()
240af33a6ddSJed Brown 
241af33a6ddSJed Brown   Level: advanced
242af33a6ddSJed Brown @*/
243bdf89e91SBarry Smith PetscErrorCode CharacteristicRegister(const char sname[],PetscErrorCode (*function)(Characteristic))
244af33a6ddSJed Brown {
245af33a6ddSJed Brown   PetscFunctionBegin;
246*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicInitializePackage());
247*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&CharacteristicList,sname,function));
248af33a6ddSJed Brown   PetscFunctionReturn(0);
249af33a6ddSJed Brown }
250af33a6ddSJed Brown 
251af33a6ddSJed Brown PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
252af33a6ddSJed Brown {
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->velocityInterp  = interp;
260af33a6ddSJed Brown   c->velocityCtx     = ctx;
261af33a6ddSJed Brown   PetscFunctionReturn(0);
262af33a6ddSJed Brown }
263af33a6ddSJed Brown 
264af33a6ddSJed Brown PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal [], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
265af33a6ddSJed Brown {
266af33a6ddSJed Brown   PetscFunctionBegin;
267af33a6ddSJed Brown   c->velocityDA          = da;
268af33a6ddSJed Brown   c->velocity            = v;
269af33a6ddSJed Brown   c->velocityOld         = vOld;
270af33a6ddSJed Brown   c->numVelocityComp     = numComponents;
271af33a6ddSJed Brown   c->velocityComp        = components;
272af33a6ddSJed Brown   c->velocityInterpLocal = interp;
273af33a6ddSJed Brown   c->velocityCtx         = ctx;
274af33a6ddSJed Brown   PetscFunctionReturn(0);
275af33a6ddSJed Brown }
276af33a6ddSJed Brown 
277af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
278af33a6ddSJed Brown {
279af33a6ddSJed Brown   PetscFunctionBegin;
280af33a6ddSJed Brown #if 0
2813c633725SBarry 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.");
282af33a6ddSJed Brown #endif
283af33a6ddSJed Brown   c->fieldDA      = da;
284af33a6ddSJed Brown   c->field        = v;
285af33a6ddSJed Brown   c->numFieldComp = numComponents;
286af33a6ddSJed Brown   c->fieldComp    = components;
287af33a6ddSJed Brown   c->fieldInterp  = interp;
288af33a6ddSJed Brown   c->fieldCtx     = ctx;
289af33a6ddSJed Brown   PetscFunctionReturn(0);
290af33a6ddSJed Brown }
291af33a6ddSJed Brown 
292af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal[], PetscInt, PetscInt[], PetscScalar [], void*), void *ctx)
293af33a6ddSJed Brown {
294af33a6ddSJed Brown   PetscFunctionBegin;
295af33a6ddSJed Brown #if 0
2963c633725SBarry 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.");
297af33a6ddSJed Brown #endif
298af33a6ddSJed Brown   c->fieldDA          = da;
299af33a6ddSJed Brown   c->field            = v;
300af33a6ddSJed Brown   c->numFieldComp     = numComponents;
301af33a6ddSJed Brown   c->fieldComp        = components;
302af33a6ddSJed Brown   c->fieldInterpLocal = interp;
303af33a6ddSJed Brown   c->fieldCtx         = ctx;
304af33a6ddSJed Brown   PetscFunctionReturn(0);
305af33a6ddSJed Brown }
306af33a6ddSJed Brown 
307af33a6ddSJed Brown PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
308af33a6ddSJed Brown {
309af33a6ddSJed Brown   CharacteristicPointDA2D Qi;
310af33a6ddSJed Brown   DM                      da = c->velocityDA;
311af33a6ddSJed Brown   Vec                     velocityLocal, velocityLocalOld;
312af33a6ddSJed Brown   Vec                     fieldLocal;
313af33a6ddSJed Brown   DMDALocalInfo           info;
314af33a6ddSJed Brown   PetscScalar             **solArray;
315af33a6ddSJed Brown   void                    *velocityArray;
316af33a6ddSJed Brown   void                    *velocityArrayOld;
317af33a6ddSJed Brown   void                    *fieldArray;
3181cc8b206SBarry Smith   PetscScalar             *interpIndices;
3191cc8b206SBarry Smith   PetscScalar             *velocityValues, *velocityValuesOld;
3201cc8b206SBarry Smith   PetscScalar             *fieldValues;
321af33a6ddSJed Brown   PetscMPIInt             rank;
322af33a6ddSJed Brown   PetscInt                dim;
323af33a6ddSJed Brown   PetscMPIInt             neighbors[9];
324af33a6ddSJed Brown   PetscInt                dof;
325af33a6ddSJed Brown   PetscInt                gx, gy;
326af33a6ddSJed Brown   PetscInt                n, is, ie, js, je, comp;
327af33a6ddSJed Brown 
328af33a6ddSJed Brown   PetscFunctionBegin;
329af33a6ddSJed Brown   c->queueSize = 0;
330*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
331*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetNeighborsRank(da, neighbors));
332*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicSetNeighbors(c, 9, neighbors));
333*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicSetUp(c));
334af33a6ddSJed Brown   /* global and local grid info */
335*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da, &dim, &gx, &gy, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
336*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da, &info));
337af33a6ddSJed Brown   is   = info.xs;          ie   = info.xs+info.xm;
338af33a6ddSJed Brown   js   = info.ys;          je   = info.ys+info.ym;
339af33a6ddSJed Brown   /* Allocation */
340*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(dim,                &interpIndices));
341*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(c->numVelocityComp, &velocityValues));
342*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(c->numVelocityComp, &velocityValuesOld));
343*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(c->numFieldComp,    &fieldValues));
344*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_Solve,NULL,NULL,NULL,NULL));
345af33a6ddSJed Brown 
346af33a6ddSJed Brown   /* -----------------------------------------------------------------------
347af33a6ddSJed Brown      PART 1, AT t-dt/2
348af33a6ddSJed Brown      -----------------------------------------------------------------------*/
349*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_QueueSetup,NULL,NULL,NULL,NULL));
350af33a6ddSJed Brown   /* GET POSITION AT HALF TIME IN THE PAST */
351af33a6ddSJed Brown   if (c->velocityInterpLocal) {
352*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(c->velocityDA, &velocityLocal));
353*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(c->velocityDA, &velocityLocalOld));
354*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
355*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
356*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
357*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
358*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(c->velocityDA, velocityLocal,    &velocityArray));
359*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
360af33a6ddSJed Brown   }
361*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n"));
362af33a6ddSJed Brown   for (Qi.j = js; Qi.j < je; Qi.j++) {
363af33a6ddSJed Brown     for (Qi.i = is; Qi.i < ie; Qi.i++) {
364af33a6ddSJed Brown       interpIndices[0] = Qi.i;
365af33a6ddSJed Brown       interpIndices[1] = Qi.j;
366*5f80ce2aSJacob Faibussowitsch       if (c->velocityInterpLocal) CHKERRQ(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
367*5f80ce2aSJacob Faibussowitsch       else CHKERRQ(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
368af33a6ddSJed Brown       Qi.x = Qi.i - velocityValues[0]*dt/2.0;
369af33a6ddSJed Brown       Qi.y = Qi.j - velocityValues[1]*dt/2.0;
370af33a6ddSJed Brown 
371af33a6ddSJed Brown       /* Determine whether the position at t - dt/2 is local */
372af33a6ddSJed Brown       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
373af33a6ddSJed Brown 
374af33a6ddSJed Brown       /* Check for Periodic boundaries and move all periodic points back onto the domain */
375*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y)));
376*5f80ce2aSJacob Faibussowitsch       CHKERRQ(CharacteristicAddPoint(c, &Qi));
377af33a6ddSJed Brown     }
378af33a6ddSJed Brown   }
379*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_QueueSetup,NULL,NULL,NULL,NULL));
380af33a6ddSJed Brown 
381*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
382*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicSendCoordinatesBegin(c));
383*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
384af33a6ddSJed Brown 
385*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,NULL,NULL,NULL,NULL));
386af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (local values) */
387*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n"));
388af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
389af33a6ddSJed Brown     Qi = c->queue[n];
390af33a6ddSJed Brown     if (c->neighbors[Qi.proc] == rank) {
391af33a6ddSJed Brown       interpIndices[0] = Qi.x;
392af33a6ddSJed Brown       interpIndices[1] = Qi.y;
393af33a6ddSJed Brown       if (c->velocityInterpLocal) {
394*5f80ce2aSJacob Faibussowitsch         CHKERRQ(c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
395*5f80ce2aSJacob Faibussowitsch         CHKERRQ(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
396af33a6ddSJed Brown       } else {
397*5f80ce2aSJacob Faibussowitsch         CHKERRQ(c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
398*5f80ce2aSJacob Faibussowitsch         CHKERRQ(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
399af33a6ddSJed Brown       }
400af33a6ddSJed Brown       Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
401af33a6ddSJed Brown       Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
402af33a6ddSJed Brown     }
403af33a6ddSJed Brown     c->queue[n] = Qi;
404af33a6ddSJed Brown   }
405*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,NULL,NULL,NULL,NULL));
406af33a6ddSJed Brown 
407*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
408*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicSendCoordinatesEnd(c));
409*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
410af33a6ddSJed Brown 
411af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (fill remote requests) */
412*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,NULL,NULL,NULL,NULL));
413*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize));
414af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
415af33a6ddSJed Brown     Qi = c->queueRemote[n];
416af33a6ddSJed Brown     interpIndices[0] = Qi.x;
417af33a6ddSJed Brown     interpIndices[1] = Qi.y;
418af33a6ddSJed Brown     if (c->velocityInterpLocal) {
419*5f80ce2aSJacob Faibussowitsch       CHKERRQ(c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx));
420*5f80ce2aSJacob Faibussowitsch       CHKERRQ(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
421af33a6ddSJed Brown     } else {
422*5f80ce2aSJacob Faibussowitsch       CHKERRQ(c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx));
423*5f80ce2aSJacob Faibussowitsch       CHKERRQ(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
424af33a6ddSJed Brown     }
425af33a6ddSJed Brown     Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
426af33a6ddSJed Brown     Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
427af33a6ddSJed Brown     c->queueRemote[n] = Qi;
428af33a6ddSJed Brown   }
429*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,NULL,NULL,NULL,NULL));
430*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
431*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicGetValuesBegin(c));
432*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicGetValuesEnd(c));
433af33a6ddSJed Brown   if (c->velocityInterpLocal) {
434*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(c->velocityDA, velocityLocal,    &velocityArray));
435*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
436*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(c->velocityDA, &velocityLocal));
437*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(c->velocityDA, &velocityLocalOld));
438af33a6ddSJed Brown   }
439*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
440af33a6ddSJed Brown 
441af33a6ddSJed Brown   /* -----------------------------------------------------------------------
442af33a6ddSJed Brown      PART 2, AT t-dt
443af33a6ddSJed Brown      -----------------------------------------------------------------------*/
444af33a6ddSJed Brown 
445af33a6ddSJed Brown   /* GET POSITION AT t_n (local values) */
446*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL));
447*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL, "Calculating position at t_{n}\n"));
448af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
449af33a6ddSJed Brown     Qi   = c->queue[n];
450af33a6ddSJed Brown     Qi.x = Qi.i - Qi.x*dt;
451af33a6ddSJed Brown     Qi.y = Qi.j - Qi.y*dt;
452af33a6ddSJed Brown 
453af33a6ddSJed Brown     /* Determine whether the position at t-dt is local */
454af33a6ddSJed Brown     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
455af33a6ddSJed Brown 
456af33a6ddSJed Brown     /* Check for Periodic boundaries and move all periodic points back onto the domain */
457*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y)));
458af33a6ddSJed Brown 
459af33a6ddSJed Brown     c->queue[n] = Qi;
460af33a6ddSJed Brown   }
461*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL));
462af33a6ddSJed Brown 
463*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
464*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicSendCoordinatesBegin(c));
465*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
466af33a6ddSJed Brown 
467af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
468*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL));
469af33a6ddSJed Brown   if (c->fieldInterpLocal) {
470*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(c->fieldDA, &fieldLocal));
471*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
472*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
473*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray));
474af33a6ddSJed Brown   }
475*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL, "Calculating local field at t_{n}\n"));
476af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
477af33a6ddSJed Brown     if (c->neighbors[c->queue[n].proc] == rank) {
478af33a6ddSJed Brown       interpIndices[0] = c->queue[n].x;
479af33a6ddSJed Brown       interpIndices[1] = c->queue[n].y;
480*5f80ce2aSJacob Faibussowitsch       if (c->fieldInterpLocal) CHKERRQ(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
481*5f80ce2aSJacob Faibussowitsch       else CHKERRQ(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
482bbd56ea5SKarl Rupp       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
483af33a6ddSJed Brown     }
484af33a6ddSJed Brown   }
485*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL));
486af33a6ddSJed Brown 
487*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
488*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicSendCoordinatesEnd(c));
489*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
490af33a6ddSJed Brown 
491af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
492*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,NULL,NULL,NULL,NULL));
493*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize));
494af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
495af33a6ddSJed Brown     interpIndices[0] = c->queueRemote[n].x;
496af33a6ddSJed Brown     interpIndices[1] = c->queueRemote[n].y;
497af33a6ddSJed Brown 
498af33a6ddSJed Brown     /* for debugging purposes */
499af33a6ddSJed Brown     if (1) { /* hacked bounds test...let's do better */
500af33a6ddSJed Brown       PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1];
501af33a6ddSJed Brown 
5023c633725SBarry Smith       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)", im, jm);
503af33a6ddSJed Brown     }
504af33a6ddSJed Brown 
505*5f80ce2aSJacob Faibussowitsch     if (c->fieldInterpLocal) CHKERRQ(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
506*5f80ce2aSJacob Faibussowitsch     else CHKERRQ(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
507bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
508af33a6ddSJed Brown   }
509*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,NULL,NULL,NULL,NULL));
510af33a6ddSJed Brown 
511*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
512*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicGetValuesBegin(c));
513*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicGetValuesEnd(c));
514af33a6ddSJed Brown   if (c->fieldInterpLocal) {
515*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray));
516*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(c->fieldDA, &fieldLocal));
517af33a6ddSJed Brown   }
518*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
519af33a6ddSJed Brown 
520af33a6ddSJed Brown   /* Return field of characteristics at t_n-1 */
521*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(CHARACTERISTIC_DAUpdate,NULL,NULL,NULL,NULL));
522*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(c->fieldDA,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
523*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(c->fieldDA, solution, &solArray));
524af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
525af33a6ddSJed Brown     Qi = c->queue[n];
526bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp];
527af33a6ddSJed Brown   }
528*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(c->fieldDA, solution, &solArray));
529*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_DAUpdate,NULL,NULL,NULL,NULL));
530*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(CHARACTERISTIC_Solve,NULL,NULL,NULL,NULL));
531af33a6ddSJed Brown 
532af33a6ddSJed Brown   /* Cleanup */
533*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(interpIndices));
534*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(velocityValues));
535*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(velocityValuesOld));
536*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(fieldValues));
537af33a6ddSJed Brown   PetscFunctionReturn(0);
538af33a6ddSJed Brown }
539af33a6ddSJed Brown 
540af33a6ddSJed Brown PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
541af33a6ddSJed Brown {
542af33a6ddSJed Brown   PetscFunctionBegin;
543af33a6ddSJed Brown   c->numNeighbors = numNeighbors;
544*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(c->neighbors));
545*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numNeighbors, &c->neighbors));
546*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(c->neighbors, neighbors, numNeighbors));
547af33a6ddSJed Brown   PetscFunctionReturn(0);
548af33a6ddSJed Brown }
549af33a6ddSJed Brown 
550af33a6ddSJed Brown PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
551af33a6ddSJed Brown {
552af33a6ddSJed Brown   PetscFunctionBegin;
5533c633725SBarry Smith   PetscCheck(c->queueSize < c->queueMax,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeded maximum queue size %d", c->queueMax);
554af33a6ddSJed Brown   c->queue[c->queueSize++] = *point;
555af33a6ddSJed Brown   PetscFunctionReturn(0);
556af33a6ddSJed Brown }
557af33a6ddSJed Brown 
558af33a6ddSJed Brown int CharacteristicSendCoordinatesBegin(Characteristic c)
559af33a6ddSJed Brown {
560af33a6ddSJed Brown   PetscMPIInt    rank, tag = 121;
561af33a6ddSJed Brown   PetscInt       i, n;
562af33a6ddSJed Brown 
563af33a6ddSJed Brown   PetscFunctionBegin;
564*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
565*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CharacteristicHeapSort(c, c->queue, c->queueSize));
566*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArrayzero(c->needCount, c->numNeighbors));
567bbd56ea5SKarl Rupp   for (i = 0;  i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
568af33a6ddSJed Brown   c->fillCount[0] = 0;
569af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
570*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1])));
571af33a6ddSJed Brown   }
572af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
573*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
574af33a6ddSJed Brown   }
575*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Waitall(c->numNeighbors-1, c->request, c->status));
576af33a6ddSJed Brown   /* Initialize the remote queue */
577af33a6ddSJed Brown   c->queueLocalMax  = c->localOffsets[0]  = 0;
578af33a6ddSJed Brown   c->queueRemoteMax = c->remoteOffsets[0] = 0;
579af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
580af33a6ddSJed Brown     c->remoteOffsets[n] = c->queueRemoteMax;
581af33a6ddSJed Brown     c->queueRemoteMax  += c->fillCount[n];
582af33a6ddSJed Brown     c->localOffsets[n]  = c->queueLocalMax;
583af33a6ddSJed Brown     c->queueLocalMax   += c->needCount[n];
584af33a6ddSJed Brown   }
585af33a6ddSJed Brown   /* HACK BEGIN */
586bbd56ea5SKarl Rupp   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
587af33a6ddSJed Brown   c->needCount[0] = 0;
588af33a6ddSJed Brown   /* HACK END */
589af33a6ddSJed Brown   if (c->queueRemoteMax) {
590*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(c->queueRemoteMax, &c->queueRemote));
5910298fd71SBarry Smith   } else c->queueRemote = NULL;
592af33a6ddSJed Brown   c->queueRemoteSize = c->queueRemoteMax;
593af33a6ddSJed Brown 
594af33a6ddSJed Brown   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
595af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
596*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]));
597*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1])));
598af33a6ddSJed Brown   }
599af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
600*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]));
601*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
602af33a6ddSJed Brown   }
603af33a6ddSJed Brown   PetscFunctionReturn(0);
604af33a6ddSJed Brown }
605af33a6ddSJed Brown 
606af33a6ddSJed Brown PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
607af33a6ddSJed Brown {
608af33a6ddSJed Brown #if 0
609af33a6ddSJed Brown   PetscMPIInt rank;
610af33a6ddSJed Brown   PetscInt    n;
611af33a6ddSJed Brown #endif
612af33a6ddSJed Brown 
613af33a6ddSJed Brown   PetscFunctionBegin;
614*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Waitall(c->numNeighbors-1, c->request, c->status));
615af33a6ddSJed Brown #if 0
616*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
617af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
6183c633725SBarry 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);
619af33a6ddSJed Brown   }
620af33a6ddSJed Brown #endif
621af33a6ddSJed Brown   PetscFunctionReturn(0);
622af33a6ddSJed Brown }
623af33a6ddSJed Brown 
624af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
625af33a6ddSJed Brown {
626af33a6ddSJed Brown   PetscMPIInt    tag = 121;
627af33a6ddSJed Brown   PetscInt       n;
628af33a6ddSJed Brown 
629af33a6ddSJed Brown   PetscFunctionBegin;
630af33a6ddSJed Brown   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
631af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
632*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1])));
633af33a6ddSJed Brown   }
634af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
635*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
636af33a6ddSJed Brown   }
637af33a6ddSJed Brown   PetscFunctionReturn(0);
638af33a6ddSJed Brown }
639af33a6ddSJed Brown 
640af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
641af33a6ddSJed Brown {
642af33a6ddSJed Brown   PetscFunctionBegin;
643*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Waitall(c->numNeighbors-1, c->request, c->status));
644af33a6ddSJed Brown   /* Free queue of requests from other procs */
645*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(c->queueRemote));
646af33a6ddSJed Brown   PetscFunctionReturn(0);
647af33a6ddSJed Brown }
648af33a6ddSJed Brown 
649af33a6ddSJed Brown /*---------------------------------------------------------------------*/
650af33a6ddSJed Brown /*
651af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
652af33a6ddSJed Brown */
6530b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
654af33a6ddSJed Brown /*---------------------------------------------------------------------*/
655af33a6ddSJed Brown {
656af33a6ddSJed Brown   CharacteristicPointDA2D temp;
6570b8f38c8SBarry Smith   PetscInt                n;
658af33a6ddSJed Brown 
6590b8f38c8SBarry Smith   PetscFunctionBegin;
660af33a6ddSJed Brown   if (0) { /* Check the order of the queue before sorting */
661*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(NULL, "Before Heap sort\n"));
662af33a6ddSJed Brown     for (n=0; n<size; n++) {
663*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(NULL,"%d %d\n",n,queue[n].proc));
664af33a6ddSJed Brown     }
665af33a6ddSJed Brown   }
666af33a6ddSJed Brown 
667af33a6ddSJed Brown   /* SORTING PHASE */
6680b8f38c8SBarry Smith   for (n = (size / 2)-1; n >= 0; n--) {
669*5f80ce2aSJacob Faibussowitsch     CHKERRQ(CharacteristicSiftDown(c, queue, n, size-1)); /* Rich had size-1 here, Matt had size*/
6700b8f38c8SBarry Smith   }
671af33a6ddSJed Brown   for (n = size-1; n >= 1; n--) {
672af33a6ddSJed Brown     temp     = queue[0];
673af33a6ddSJed Brown     queue[0] = queue[n];
674af33a6ddSJed Brown     queue[n] = temp;
675*5f80ce2aSJacob Faibussowitsch     CHKERRQ(CharacteristicSiftDown(c, queue, 0, n-1));
676af33a6ddSJed Brown   }
677af33a6ddSJed Brown   if (0) { /* Check the order of the queue after sorting */
678*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(NULL, "Avter  Heap sort\n"));
679af33a6ddSJed Brown     for (n=0; n<size; n++) {
680*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(NULL,"%d %d\n",n,queue[n].proc));
681af33a6ddSJed Brown     }
682af33a6ddSJed Brown   }
6830b8f38c8SBarry Smith   PetscFunctionReturn(0);
684af33a6ddSJed Brown }
685af33a6ddSJed Brown 
686af33a6ddSJed Brown /*---------------------------------------------------------------------*/
687af33a6ddSJed Brown /*
688af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
689af33a6ddSJed Brown */
6900b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
691af33a6ddSJed Brown /*---------------------------------------------------------------------*/
692af33a6ddSJed Brown {
693af33a6ddSJed Brown   PetscBool               done = PETSC_FALSE;
694af33a6ddSJed Brown   PetscInt                maxChild;
695af33a6ddSJed Brown   CharacteristicPointDA2D temp;
696af33a6ddSJed Brown 
697af33a6ddSJed Brown   PetscFunctionBegin;
698af33a6ddSJed Brown   while ((root*2 <= bottom) && (!done)) {
699af33a6ddSJed Brown     if (root*2 == bottom) maxChild = root * 2;
700af33a6ddSJed Brown     else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2;
701af33a6ddSJed Brown     else maxChild = root * 2 + 1;
702af33a6ddSJed Brown 
703af33a6ddSJed Brown     if (queue[root].proc < queue[maxChild].proc) {
704af33a6ddSJed Brown       temp            = queue[root];
705af33a6ddSJed Brown       queue[root]     = queue[maxChild];
706af33a6ddSJed Brown       queue[maxChild] = temp;
707af33a6ddSJed Brown       root            = maxChild;
708db4deed7SKarl Rupp     } else done = PETSC_TRUE;
709af33a6ddSJed Brown   }
710af33a6ddSJed Brown   PetscFunctionReturn(0);
711af33a6ddSJed Brown }
712af33a6ddSJed Brown 
713af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
714af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
715af33a6ddSJed Brown {
716bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx, by;
717af33a6ddSJed Brown   PetscBool        IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
718af33a6ddSJed Brown   MPI_Comm         comm;
719af33a6ddSJed Brown   PetscMPIInt      rank;
720af33a6ddSJed Brown   PetscInt         **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ;
721af33a6ddSJed Brown 
722af33a6ddSJed Brown   PetscFunctionBegin;
723*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) da, &comm));
724*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
725*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da, NULL, NULL, NULL, NULL, &PI,&PJ, NULL, NULL, NULL, &bx, &by,NULL, NULL));
726af33a6ddSJed Brown 
727bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
728bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
729af33a6ddSJed Brown 
730af33a6ddSJed Brown   neighbors[0] = rank;
731af33a6ddSJed Brown   rank = 0;
732*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(PJ,&procs));
733af33a6ddSJed Brown   for (pj=0; pj<PJ; pj++) {
734*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(PI,&(procs[pj])));
735af33a6ddSJed Brown     for (pi=0; pi<PI; pi++) {
736af33a6ddSJed Brown       procs[pj][pi] = rank;
737af33a6ddSJed Brown       rank++;
738af33a6ddSJed Brown     }
739af33a6ddSJed Brown   }
740af33a6ddSJed Brown 
741af33a6ddSJed Brown   pi  = neighbors[0] % PI;
742af33a6ddSJed Brown   pj  = neighbors[0] / PI;
743af33a6ddSJed Brown   pim = pi-1;  if (pim<0) pim=PI-1;
744af33a6ddSJed Brown   pip = (pi+1)%PI;
745af33a6ddSJed Brown   pjm = pj-1;  if (pjm<0) pjm=PJ-1;
746af33a6ddSJed Brown   pjp = (pj+1)%PJ;
747af33a6ddSJed Brown 
748af33a6ddSJed Brown   neighbors[1] = procs[pj] [pim];
749af33a6ddSJed Brown   neighbors[2] = procs[pjp][pim];
750af33a6ddSJed Brown   neighbors[3] = procs[pjp][pi];
751af33a6ddSJed Brown   neighbors[4] = procs[pjp][pip];
752af33a6ddSJed Brown   neighbors[5] = procs[pj] [pip];
753af33a6ddSJed Brown   neighbors[6] = procs[pjm][pip];
754af33a6ddSJed Brown   neighbors[7] = procs[pjm][pi];
755af33a6ddSJed Brown   neighbors[8] = procs[pjm][pim];
756af33a6ddSJed Brown 
757af33a6ddSJed Brown   if (!IPeriodic) {
758af33a6ddSJed Brown     if (pi==0)    neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0];
759af33a6ddSJed Brown     if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0];
760af33a6ddSJed Brown   }
761af33a6ddSJed Brown 
762af33a6ddSJed Brown   if (!JPeriodic) {
763af33a6ddSJed Brown     if (pj==0)    neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0];
764af33a6ddSJed Brown     if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0];
765af33a6ddSJed Brown   }
766af33a6ddSJed Brown 
767af33a6ddSJed Brown   for (pj = 0; pj < PJ; pj++) {
768*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(procs[pj]));
769af33a6ddSJed Brown   }
770*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(procs));
771af33a6ddSJed Brown   PetscFunctionReturn(0);
772af33a6ddSJed Brown }
773af33a6ddSJed Brown 
774af33a6ddSJed Brown /*
775af33a6ddSJed Brown   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
776af33a6ddSJed Brown     2 | 3 | 4
777af33a6ddSJed Brown     __|___|__
778af33a6ddSJed Brown     1 | 0 | 5
779af33a6ddSJed Brown     __|___|__
780af33a6ddSJed Brown     8 | 7 | 6
781af33a6ddSJed Brown       |   |
782af33a6ddSJed Brown */
7831cc8b206SBarry Smith PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr)
784af33a6ddSJed Brown {
785af33a6ddSJed Brown   DMDALocalInfo  info;
7861cc8b206SBarry Smith   PetscReal      is,ie,js,je;
787af33a6ddSJed Brown 
788*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da, &info));
7891cc8b206SBarry Smith   is   = (PetscReal) info.xs - 0.5; ie = (PetscReal) info.xs + info.xm - 0.5;
7901cc8b206SBarry Smith   js   = (PetscReal) info.ys - 0.5; je = (PetscReal) info.ys + info.ym - 0.5;
791af33a6ddSJed Brown 
792af33a6ddSJed Brown   if (ir >= is && ir <= ie) { /* center column */
793bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 0;
794bbd56ea5SKarl Rupp     else if (jr < js)         return 7;
795bbd56ea5SKarl Rupp     else                      return 3;
796af33a6ddSJed Brown   } else if (ir < is) {     /* left column */
797bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 1;
798bbd56ea5SKarl Rupp     else if (jr < js)         return 8;
799bbd56ea5SKarl Rupp     else                      return 2;
800af33a6ddSJed Brown   } else {                  /* right column */
801bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 5;
802bbd56ea5SKarl Rupp     else if (jr < js)         return 6;
803bbd56ea5SKarl Rupp     else                      return 4;
804af33a6ddSJed Brown   }
805af33a6ddSJed Brown }
806