xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision 1e25c2744a56dc779d05d2a7975a7f2272f1aa7b)
1af33a6ddSJed Brown 
2b45d2f2cSJed Brown #include <petsc-private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/
3*1e25c274SJed Brown #include <petscdmda.h>
4af33a6ddSJed Brown 
5af33a6ddSJed Brown PetscClassId  CHARACTERISTIC_CLASSID;
6af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate;
7af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange;
8af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange;
9af33a6ddSJed Brown /*
10af33a6ddSJed Brown    Contains the list of registered characteristic routines
11af33a6ddSJed Brown */
120298fd71SBarry Smith PetscFunctionList CharacteristicList              = NULL;
13140e18c1SBarry Smith PetscBool         CharacteristicRegisterAllCalled = PETSC_FALSE;
14af33a6ddSJed Brown 
15af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt []);
16af33a6ddSJed Brown PetscInt       DMDAGetNeighborRelative(DM, PassiveReal, PassiveReal);
17af33a6ddSJed Brown PetscErrorCode DMDAMapToPeriodicDomain(DM, PetscScalar []);
18af33a6ddSJed Brown 
190b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt);
200b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt);
21af33a6ddSJed Brown 
22af33a6ddSJed Brown #undef __FUNCT__
23af33a6ddSJed Brown #define __FUNCT__ "CharacteristicView"
24af33a6ddSJed Brown PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
25af33a6ddSJed Brown {
26af33a6ddSJed Brown   PetscBool      iascii;
27af33a6ddSJed Brown   PetscErrorCode ierr;
28af33a6ddSJed Brown 
29af33a6ddSJed Brown   PetscFunctionBegin;
30af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
31af33a6ddSJed Brown   if (!viewer) {
32ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr);
33af33a6ddSJed Brown   }
34af33a6ddSJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
35af33a6ddSJed Brown   PetscCheckSameComm(c, 1, viewer, 2);
36af33a6ddSJed Brown 
37251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
38af33a6ddSJed Brown   if (!iascii) {
39af33a6ddSJed Brown     if (c->ops->view) {
40af33a6ddSJed Brown       ierr = (*c->ops->view)(c, viewer);CHKERRQ(ierr);
41af33a6ddSJed Brown     }
42af33a6ddSJed Brown   }
43af33a6ddSJed Brown   PetscFunctionReturn(0);
44af33a6ddSJed Brown }
45af33a6ddSJed Brown 
46af33a6ddSJed Brown #undef __FUNCT__
47af33a6ddSJed Brown #define __FUNCT__ "CharacteristicDestroy"
486bf464f9SBarry Smith PetscErrorCode CharacteristicDestroy(Characteristic *c)
49af33a6ddSJed Brown {
50af33a6ddSJed Brown   PetscErrorCode ierr;
51af33a6ddSJed Brown 
52af33a6ddSJed Brown   PetscFunctionBegin;
536bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
546bf464f9SBarry Smith   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
556bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0);
56af33a6ddSJed Brown 
576bf464f9SBarry Smith   if ((*c)->ops->destroy) {
586bf464f9SBarry Smith     ierr = (*(*c)->ops->destroy)((*c));CHKERRQ(ierr);
59af33a6ddSJed Brown   }
606bf464f9SBarry Smith   ierr = MPI_Type_free(&(*c)->itemType);CHKERRQ(ierr);
616bf464f9SBarry Smith   ierr = PetscFree((*c)->queue);CHKERRQ(ierr);
626bf464f9SBarry Smith   ierr = PetscFree((*c)->queueLocal);CHKERRQ(ierr);
636bf464f9SBarry Smith   ierr = PetscFree((*c)->queueRemote);CHKERRQ(ierr);
646bf464f9SBarry Smith   ierr = PetscFree((*c)->neighbors);CHKERRQ(ierr);
656bf464f9SBarry Smith   ierr = PetscFree((*c)->needCount);CHKERRQ(ierr);
666bf464f9SBarry Smith   ierr = PetscFree((*c)->localOffsets);CHKERRQ(ierr);
676bf464f9SBarry Smith   ierr = PetscFree((*c)->fillCount);CHKERRQ(ierr);
686bf464f9SBarry Smith   ierr = PetscFree((*c)->remoteOffsets);CHKERRQ(ierr);
696bf464f9SBarry Smith   ierr = PetscFree((*c)->request);CHKERRQ(ierr);
706bf464f9SBarry Smith   ierr = PetscFree((*c)->status);CHKERRQ(ierr);
71af33a6ddSJed Brown   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
72af33a6ddSJed Brown   PetscFunctionReturn(0);
73af33a6ddSJed Brown }
74af33a6ddSJed Brown 
75af33a6ddSJed Brown #undef __FUNCT__
76af33a6ddSJed Brown #define __FUNCT__ "CharacteristicCreate"
77af33a6ddSJed Brown PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
78af33a6ddSJed Brown {
79af33a6ddSJed Brown   Characteristic newC;
80af33a6ddSJed Brown   PetscErrorCode ierr;
81af33a6ddSJed Brown 
82af33a6ddSJed Brown   PetscFunctionBegin;
83af33a6ddSJed Brown   PetscValidPointer(c, 2);
840298fd71SBarry Smith   *c = NULL;
85519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
860298fd71SBarry Smith   ierr = CharacteristicInitializePackage(NULL);CHKERRQ(ierr);
87af33a6ddSJed Brown #endif
88af33a6ddSJed Brown 
8967c2884eSBarry Smith   ierr = PetscHeaderCreate(newC, _p_Characteristic, struct _CharacteristicOps, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "SemiLagrange", comm, CharacteristicDestroy, CharacteristicView);CHKERRQ(ierr);
90af33a6ddSJed Brown   ierr = PetscLogObjectCreate(newC);CHKERRQ(ierr);
91af33a6ddSJed Brown   *c   = newC;
92af33a6ddSJed Brown 
93af33a6ddSJed Brown   newC->structured          = PETSC_TRUE;
94af33a6ddSJed Brown   newC->numIds              = 0;
950298fd71SBarry Smith   newC->velocityDA          = NULL;
960298fd71SBarry Smith   newC->velocity            = NULL;
970298fd71SBarry Smith   newC->velocityOld         = NULL;
98af33a6ddSJed Brown   newC->numVelocityComp     = 0;
990298fd71SBarry Smith   newC->velocityComp        = NULL;
1000298fd71SBarry Smith   newC->velocityInterp      = NULL;
1010298fd71SBarry Smith   newC->velocityInterpLocal = NULL;
1020298fd71SBarry Smith   newC->velocityCtx         = NULL;
1030298fd71SBarry Smith   newC->fieldDA             = NULL;
1040298fd71SBarry Smith   newC->field               = NULL;
105af33a6ddSJed Brown   newC->numFieldComp        = 0;
1060298fd71SBarry Smith   newC->fieldComp           = NULL;
1070298fd71SBarry Smith   newC->fieldInterp         = NULL;
1080298fd71SBarry Smith   newC->fieldInterpLocal    = NULL;
1090298fd71SBarry Smith   newC->fieldCtx            = NULL;
1100298fd71SBarry Smith   newC->itemType            = 0;
1110298fd71SBarry Smith   newC->queue               = NULL;
112af33a6ddSJed Brown   newC->queueSize           = 0;
113af33a6ddSJed Brown   newC->queueMax            = 0;
1140298fd71SBarry Smith   newC->queueLocal          = NULL;
115af33a6ddSJed Brown   newC->queueLocalSize      = 0;
116af33a6ddSJed Brown   newC->queueLocalMax       = 0;
1170298fd71SBarry Smith   newC->queueRemote         = NULL;
118af33a6ddSJed Brown   newC->queueRemoteSize     = 0;
119af33a6ddSJed Brown   newC->queueRemoteMax      = 0;
120af33a6ddSJed Brown   newC->numNeighbors        = 0;
1210298fd71SBarry Smith   newC->neighbors           = NULL;
1220298fd71SBarry Smith   newC->needCount           = NULL;
1230298fd71SBarry Smith   newC->localOffsets        = NULL;
1240298fd71SBarry Smith   newC->fillCount           = NULL;
1250298fd71SBarry Smith   newC->remoteOffsets       = NULL;
1260298fd71SBarry Smith   newC->request             = NULL;
1270298fd71SBarry Smith   newC->status              = NULL;
128af33a6ddSJed Brown   PetscFunctionReturn(0);
129af33a6ddSJed Brown }
130af33a6ddSJed Brown 
131af33a6ddSJed Brown #undef __FUNCT__
132af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetType"
133af33a6ddSJed Brown /*@C
134af33a6ddSJed Brown    CharacteristicSetType - Builds Characteristic for a particular solver.
135af33a6ddSJed Brown 
136af33a6ddSJed Brown    Logically Collective on Characteristic
137af33a6ddSJed Brown 
138af33a6ddSJed Brown    Input Parameters:
139af33a6ddSJed Brown +  c    - the method of characteristics context
140af33a6ddSJed Brown -  type - a known method
141af33a6ddSJed Brown 
142af33a6ddSJed Brown    Options Database Key:
143af33a6ddSJed Brown .  -characteristic_type <method> - Sets the method; use -help for a list
144af33a6ddSJed Brown     of available methods
145af33a6ddSJed Brown 
146af33a6ddSJed Brown    Notes:
14730a76a96SBarry Smith    See "include/petsccharacteristic.h" for available methods
148af33a6ddSJed Brown 
149af33a6ddSJed Brown   Normally, it is best to use the CharacteristicSetFromOptions() command and
150af33a6ddSJed Brown   then set the Characteristic type from the options database rather than by using
151af33a6ddSJed Brown   this routine.  Using the options database provides the user with
152af33a6ddSJed Brown   maximum flexibility in evaluating the many different Krylov methods.
153af33a6ddSJed Brown   The CharacteristicSetType() routine is provided for those situations where it
154af33a6ddSJed Brown   is necessary to set the iterative solver independently of the command
155af33a6ddSJed Brown   line or options database.  This might be the case, for example, when
156af33a6ddSJed Brown   the choice of iterative solver changes during the execution of the
157af33a6ddSJed Brown   program, and the user's application is taking responsibility for
158af33a6ddSJed Brown   choosing the appropriate method.  In other words, this routine is
159af33a6ddSJed Brown   not for beginners.
160af33a6ddSJed Brown 
161af33a6ddSJed Brown   Level: intermediate
162af33a6ddSJed Brown 
163af33a6ddSJed Brown .keywords: Characteristic, set, method
164af33a6ddSJed Brown 
165af33a6ddSJed Brown .seealso: CharacteristicType
166af33a6ddSJed Brown 
167af33a6ddSJed Brown @*/
16819fd82e9SBarry Smith PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
169af33a6ddSJed Brown {
170af33a6ddSJed Brown   PetscErrorCode ierr, (*r)(Characteristic);
171af33a6ddSJed Brown   PetscBool      match;
172af33a6ddSJed Brown 
173af33a6ddSJed Brown   PetscFunctionBegin;
174af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
175af33a6ddSJed Brown   PetscValidCharPointer(type, 2);
176af33a6ddSJed Brown 
177251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) c, type, &match);CHKERRQ(ierr);
178af33a6ddSJed Brown   if (match) PetscFunctionReturn(0);
179af33a6ddSJed Brown 
180af33a6ddSJed Brown   if (c->data) {
181af33a6ddSJed Brown     /* destroy the old private Characteristic context */
182af33a6ddSJed Brown     ierr = (*c->ops->destroy)(c);CHKERRQ(ierr);
1830298fd71SBarry Smith     c->ops->destroy = NULL;
184af33a6ddSJed Brown     c->data         = 0;
185af33a6ddSJed Brown   }
186af33a6ddSJed Brown 
187ce94432eSBarry Smith   ierr =  PetscFunctionListFind(PetscObjectComm((PetscObject)c),CharacteristicList, type,PETSC_TRUE, (void (**)(void)) &r);CHKERRQ(ierr);
188af33a6ddSJed Brown   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
189af33a6ddSJed Brown   c->setupcalled = 0;
190af33a6ddSJed Brown   ierr = (*r)(c);CHKERRQ(ierr);
191af33a6ddSJed Brown   ierr = PetscObjectChangeTypeName((PetscObject) c, type);CHKERRQ(ierr);
192af33a6ddSJed Brown   PetscFunctionReturn(0);
193af33a6ddSJed Brown }
194af33a6ddSJed Brown 
195af33a6ddSJed Brown #undef __FUNCT__
196af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetUp"
197af33a6ddSJed Brown /*@
198af33a6ddSJed Brown    CharacteristicSetUp - Sets up the internal data structures for the
199af33a6ddSJed Brown    later use of an iterative solver.
200af33a6ddSJed Brown 
201af33a6ddSJed Brown    Collective on Characteristic
202af33a6ddSJed Brown 
203af33a6ddSJed Brown    Input Parameter:
204af33a6ddSJed Brown .  ksp   - iterative context obtained from CharacteristicCreate()
205af33a6ddSJed Brown 
206af33a6ddSJed Brown    Level: developer
207af33a6ddSJed Brown 
208af33a6ddSJed Brown .keywords: Characteristic, setup
209af33a6ddSJed Brown 
210af33a6ddSJed Brown .seealso: CharacteristicCreate(), CharacteristicSolve(), CharacteristicDestroy()
211af33a6ddSJed Brown @*/
212af33a6ddSJed Brown PetscErrorCode CharacteristicSetUp(Characteristic c)
213af33a6ddSJed Brown {
214af33a6ddSJed Brown   PetscErrorCode ierr;
215af33a6ddSJed Brown 
216af33a6ddSJed Brown   PetscFunctionBegin;
217af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
218af33a6ddSJed Brown 
219af33a6ddSJed Brown   if (!((PetscObject)c)->type_name) {
220af33a6ddSJed Brown     ierr = CharacteristicSetType(c, CHARACTERISTICDA);CHKERRQ(ierr);
221af33a6ddSJed Brown   }
222af33a6ddSJed Brown 
223af33a6ddSJed Brown   if (c->setupcalled == 2) PetscFunctionReturn(0);
224af33a6ddSJed Brown 
225af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr);
226af33a6ddSJed Brown   if (!c->setupcalled) {
227af33a6ddSJed Brown     ierr = (*c->ops->setup)(c);CHKERRQ(ierr);
228af33a6ddSJed Brown   }
229af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr);
230af33a6ddSJed Brown   c->setupcalled = 2;
231af33a6ddSJed Brown   PetscFunctionReturn(0);
232af33a6ddSJed Brown }
233af33a6ddSJed Brown 
234af33a6ddSJed Brown #undef __FUNCT__
235af33a6ddSJed Brown #define __FUNCT__ "CharacteristicRegister"
236af33a6ddSJed Brown /*@C
237af33a6ddSJed Brown   CharacteristicRegister - See CharacteristicRegisterDynamic()
238af33a6ddSJed Brown 
239af33a6ddSJed Brown   Level: advanced
240af33a6ddSJed Brown @*/
241af33a6ddSJed Brown PetscErrorCode CharacteristicRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Characteristic))
242af33a6ddSJed Brown {
243af33a6ddSJed Brown   PetscErrorCode ierr;
244af33a6ddSJed Brown   char           fullname[PETSC_MAX_PATH_LEN];
245af33a6ddSJed Brown 
246af33a6ddSJed Brown   PetscFunctionBegin;
247140e18c1SBarry Smith   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
248140e18c1SBarry Smith   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&CharacteristicList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
249af33a6ddSJed Brown   PetscFunctionReturn(0);
250af33a6ddSJed Brown }
251af33a6ddSJed Brown 
252af33a6ddSJed Brown #undef __FUNCT__
253af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetVelocityInterpolation"
254af33a6ddSJed 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)
255af33a6ddSJed Brown {
256af33a6ddSJed Brown   PetscFunctionBegin;
257af33a6ddSJed Brown   c->velocityDA      = da;
258af33a6ddSJed Brown   c->velocity        = v;
259af33a6ddSJed Brown   c->velocityOld     = vOld;
260af33a6ddSJed Brown   c->numVelocityComp = numComponents;
261af33a6ddSJed Brown   c->velocityComp    = components;
262af33a6ddSJed Brown   c->velocityInterp  = interp;
263af33a6ddSJed Brown   c->velocityCtx     = ctx;
264af33a6ddSJed Brown   PetscFunctionReturn(0);
265af33a6ddSJed Brown }
266af33a6ddSJed Brown 
267af33a6ddSJed Brown #undef __FUNCT__
268af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetVelocityInterpolationLocal"
269af33a6ddSJed 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)
270af33a6ddSJed Brown {
271af33a6ddSJed Brown   PetscFunctionBegin;
272af33a6ddSJed Brown   c->velocityDA          = da;
273af33a6ddSJed Brown   c->velocity            = v;
274af33a6ddSJed Brown   c->velocityOld         = vOld;
275af33a6ddSJed Brown   c->numVelocityComp     = numComponents;
276af33a6ddSJed Brown   c->velocityComp        = components;
277af33a6ddSJed Brown   c->velocityInterpLocal = interp;
278af33a6ddSJed Brown   c->velocityCtx         = ctx;
279af33a6ddSJed Brown   PetscFunctionReturn(0);
280af33a6ddSJed Brown }
281af33a6ddSJed Brown 
282af33a6ddSJed Brown #undef __FUNCT__
283af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetFieldInterpolation"
284af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
285af33a6ddSJed Brown {
286af33a6ddSJed Brown   PetscFunctionBegin;
287af33a6ddSJed Brown #if 0
288af33a6ddSJed Brown   if (numComponents > 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
289af33a6ddSJed Brown #endif
290af33a6ddSJed Brown   c->fieldDA      = da;
291af33a6ddSJed Brown   c->field        = v;
292af33a6ddSJed Brown   c->numFieldComp = numComponents;
293af33a6ddSJed Brown   c->fieldComp    = components;
294af33a6ddSJed Brown   c->fieldInterp  = interp;
295af33a6ddSJed Brown   c->fieldCtx     = ctx;
296af33a6ddSJed Brown   PetscFunctionReturn(0);
297af33a6ddSJed Brown }
298af33a6ddSJed Brown 
299af33a6ddSJed Brown #undef __FUNCT__
300af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetFieldInterpolationLocal"
301af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal[], PetscInt, PetscInt[], PetscScalar [], void*), void *ctx)
302af33a6ddSJed Brown {
303af33a6ddSJed Brown   PetscFunctionBegin;
304af33a6ddSJed Brown #if 0
305af33a6ddSJed Brown   if (numComponents > 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
306af33a6ddSJed Brown #endif
307af33a6ddSJed Brown   c->fieldDA          = da;
308af33a6ddSJed Brown   c->field            = v;
309af33a6ddSJed Brown   c->numFieldComp     = numComponents;
310af33a6ddSJed Brown   c->fieldComp        = components;
311af33a6ddSJed Brown   c->fieldInterpLocal = interp;
312af33a6ddSJed Brown   c->fieldCtx         = ctx;
313af33a6ddSJed Brown   PetscFunctionReturn(0);
314af33a6ddSJed Brown }
315af33a6ddSJed Brown 
316af33a6ddSJed Brown #undef __FUNCT__
317af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSolve"
318af33a6ddSJed Brown PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
319af33a6ddSJed Brown {
320af33a6ddSJed Brown   CharacteristicPointDA2D Qi;
321af33a6ddSJed Brown   DM                      da = c->velocityDA;
322af33a6ddSJed Brown   Vec                     velocityLocal, velocityLocalOld;
323af33a6ddSJed Brown   Vec                     fieldLocal;
324af33a6ddSJed Brown   DMDALocalInfo           info;
325af33a6ddSJed Brown   PetscScalar             **solArray;
326af33a6ddSJed Brown   void                    *velocityArray;
327af33a6ddSJed Brown   void                    *velocityArrayOld;
328af33a6ddSJed Brown   void                    *fieldArray;
329af33a6ddSJed Brown   PassiveScalar           *interpIndices;
330af33a6ddSJed Brown   PassiveScalar           *velocityValues, *velocityValuesOld;
331af33a6ddSJed Brown   PassiveScalar           *fieldValues;
332af33a6ddSJed Brown   PetscMPIInt             rank;
333af33a6ddSJed Brown   PetscInt                dim;
334af33a6ddSJed Brown   PetscMPIInt             neighbors[9];
335af33a6ddSJed Brown   PetscInt                dof;
336af33a6ddSJed Brown   PetscInt                gx, gy;
337af33a6ddSJed Brown   PetscInt                n, is, ie, js, je, comp;
338af33a6ddSJed Brown   PetscErrorCode          ierr;
339af33a6ddSJed Brown 
340af33a6ddSJed Brown   PetscFunctionBegin;
341af33a6ddSJed Brown   c->queueSize = 0;
342ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
343af33a6ddSJed Brown   ierr = DMDAGetNeighborsRank(da, neighbors);CHKERRQ(ierr);
344af33a6ddSJed Brown   ierr = CharacteristicSetNeighbors(c, 9, neighbors);CHKERRQ(ierr);
345af33a6ddSJed Brown   ierr = CharacteristicSetUp(c);CHKERRQ(ierr);
346af33a6ddSJed Brown   /* global and local grid info */
347af33a6ddSJed Brown   ierr = DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);CHKERRQ(ierr);
348af33a6ddSJed Brown   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
349af33a6ddSJed Brown   is   = info.xs;          ie   = info.xs+info.xm;
350af33a6ddSJed Brown   js   = info.ys;          je   = info.ys+info.ym;
351af33a6ddSJed Brown   /* Allocation */
352af33a6ddSJed Brown   ierr = PetscMalloc(dim*sizeof(PetscScalar),                &interpIndices);CHKERRQ(ierr);
353af33a6ddSJed Brown   ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValues);CHKERRQ(ierr);
354af33a6ddSJed Brown   ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValuesOld);CHKERRQ(ierr);
355af33a6ddSJed Brown   ierr = PetscMalloc(c->numFieldComp*sizeof(PetscScalar),    &fieldValues);CHKERRQ(ierr);
356af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
357af33a6ddSJed Brown 
358af33a6ddSJed Brown   /* -----------------------------------------------------------------------
359af33a6ddSJed Brown      PART 1, AT t-dt/2
360af33a6ddSJed Brown      -----------------------------------------------------------------------*/
361af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
362af33a6ddSJed Brown   /* GET POSITION AT HALF TIME IN THE PAST */
363af33a6ddSJed Brown   if (c->velocityInterpLocal) {
364af33a6ddSJed Brown     ierr = DMGetLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
365af33a6ddSJed Brown     ierr = DMGetLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
366af33a6ddSJed Brown     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
367af33a6ddSJed Brown     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
368af33a6ddSJed Brown     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
369af33a6ddSJed Brown     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
370af33a6ddSJed Brown     ierr = DMDAVecGetArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
371af33a6ddSJed Brown     ierr = DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
372af33a6ddSJed Brown   }
3730298fd71SBarry Smith   ierr = PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n");CHKERRQ(ierr);
374af33a6ddSJed Brown   for (Qi.j = js; Qi.j < je; Qi.j++) {
375af33a6ddSJed Brown     for (Qi.i = is; Qi.i < ie; Qi.i++) {
376af33a6ddSJed Brown       interpIndices[0] = Qi.i;
377af33a6ddSJed Brown       interpIndices[1] = Qi.j;
378bbd56ea5SKarl Rupp       if (c->velocityInterpLocal) c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
379bbd56ea5SKarl Rupp       else c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
380af33a6ddSJed Brown       Qi.x = Qi.i - velocityValues[0]*dt/2.0;
381af33a6ddSJed Brown       Qi.y = Qi.j - velocityValues[1]*dt/2.0;
382af33a6ddSJed Brown 
383af33a6ddSJed Brown       /* Determine whether the position at t - dt/2 is local */
384af33a6ddSJed Brown       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
385af33a6ddSJed Brown 
386af33a6ddSJed Brown       /* Check for Periodic boundaries and move all periodic points back onto the domain */
387af33a6ddSJed Brown       ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
388af33a6ddSJed Brown       ierr = CharacteristicAddPoint(c, &Qi);CHKERRQ(ierr);
389af33a6ddSJed Brown     }
390af33a6ddSJed Brown   }
391af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
392af33a6ddSJed Brown 
393af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
394af33a6ddSJed Brown   ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr);
395af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
396af33a6ddSJed Brown 
397af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
398af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (local values) */
3990298fd71SBarry Smith   ierr = PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n");CHKERRQ(ierr);
400af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
401af33a6ddSJed Brown     Qi = c->queue[n];
402af33a6ddSJed Brown     if (c->neighbors[Qi.proc] == rank) {
403af33a6ddSJed Brown       interpIndices[0] = Qi.x;
404af33a6ddSJed Brown       interpIndices[1] = Qi.y;
405af33a6ddSJed Brown       if (c->velocityInterpLocal) {
406af33a6ddSJed Brown         c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
407af33a6ddSJed Brown         c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
408af33a6ddSJed Brown       } else {
409af33a6ddSJed Brown         c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
410af33a6ddSJed Brown         c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
411af33a6ddSJed Brown       }
412af33a6ddSJed Brown       Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
413af33a6ddSJed Brown       Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
414af33a6ddSJed Brown     }
415af33a6ddSJed Brown     c->queue[n] = Qi;
416af33a6ddSJed Brown   }
417af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
418af33a6ddSJed Brown 
419af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
420af33a6ddSJed Brown   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
421af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
422af33a6ddSJed Brown 
423af33a6ddSJed Brown 
424af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (fill remote requests) */
425af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
4260298fd71SBarry Smith   ierr = PetscInfo1(NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);CHKERRQ(ierr);
427af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
428af33a6ddSJed Brown     Qi = c->queueRemote[n];
429af33a6ddSJed Brown     interpIndices[0] = Qi.x;
430af33a6ddSJed Brown     interpIndices[1] = Qi.y;
431af33a6ddSJed Brown     if (c->velocityInterpLocal) {
432af33a6ddSJed Brown       c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
433af33a6ddSJed Brown       c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
434af33a6ddSJed Brown     } else {
435af33a6ddSJed Brown       c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
436af33a6ddSJed Brown       c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
437af33a6ddSJed Brown     }
438af33a6ddSJed Brown     Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
439af33a6ddSJed Brown     Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
440af33a6ddSJed Brown     c->queueRemote[n] = Qi;
441af33a6ddSJed Brown   }
442af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
443af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
444af33a6ddSJed Brown   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
445af33a6ddSJed Brown   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
446af33a6ddSJed Brown   if (c->velocityInterpLocal) {
447af33a6ddSJed Brown     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
448af33a6ddSJed Brown     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
449af33a6ddSJed Brown     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
450af33a6ddSJed Brown     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
451af33a6ddSJed Brown   }
452af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
453af33a6ddSJed Brown 
454af33a6ddSJed Brown   /* -----------------------------------------------------------------------
455af33a6ddSJed Brown      PART 2, AT t-dt
456af33a6ddSJed Brown      -----------------------------------------------------------------------*/
457af33a6ddSJed Brown 
458af33a6ddSJed Brown   /* GET POSITION AT t_n (local values) */
459af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
4600298fd71SBarry Smith   ierr = PetscInfo(NULL, "Calculating position at t_{n}\n");CHKERRQ(ierr);
461af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
462af33a6ddSJed Brown     Qi   = c->queue[n];
463af33a6ddSJed Brown     Qi.x = Qi.i - Qi.x*dt;
464af33a6ddSJed Brown     Qi.y = Qi.j - Qi.y*dt;
465af33a6ddSJed Brown 
466af33a6ddSJed Brown     /* Determine whether the position at t-dt is local */
467af33a6ddSJed Brown     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
468af33a6ddSJed Brown 
469af33a6ddSJed Brown     /* Check for Periodic boundaries and move all periodic points back onto the domain */
470af33a6ddSJed Brown     ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
471af33a6ddSJed Brown 
472af33a6ddSJed Brown     c->queue[n] = Qi;
473af33a6ddSJed Brown   }
474af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
475af33a6ddSJed Brown 
476af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
477af33a6ddSJed Brown   ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr);
478af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
479af33a6ddSJed Brown 
480af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
481af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
482af33a6ddSJed Brown   if (c->fieldInterpLocal) {
483af33a6ddSJed Brown     ierr = DMGetLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
484af33a6ddSJed Brown     ierr = DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
485af33a6ddSJed Brown     ierr = DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
486af33a6ddSJed Brown     ierr = DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
487af33a6ddSJed Brown   }
4880298fd71SBarry Smith   ierr = PetscInfo(NULL, "Calculating local field at t_{n}\n");CHKERRQ(ierr);
489af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
490af33a6ddSJed Brown     if (c->neighbors[c->queue[n].proc] == rank) {
491af33a6ddSJed Brown       interpIndices[0] = c->queue[n].x;
492af33a6ddSJed Brown       interpIndices[1] = c->queue[n].y;
493bbd56ea5SKarl Rupp       if (c->fieldInterpLocal) c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
494bbd56ea5SKarl Rupp       else c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
495bbd56ea5SKarl Rupp       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
496af33a6ddSJed Brown     }
497af33a6ddSJed Brown   }
498af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
499af33a6ddSJed Brown 
500af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
501af33a6ddSJed Brown   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
502af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
503af33a6ddSJed Brown 
504af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
505af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
5060298fd71SBarry Smith   ierr = PetscInfo1(NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);CHKERRQ(ierr);
507af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
508af33a6ddSJed Brown     interpIndices[0] = c->queueRemote[n].x;
509af33a6ddSJed Brown     interpIndices[1] = c->queueRemote[n].y;
510af33a6ddSJed Brown 
511af33a6ddSJed Brown     /* for debugging purposes */
512af33a6ddSJed Brown     if (1) { /* hacked bounds test...let's do better */
513af33a6ddSJed Brown       PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1];
514af33a6ddSJed Brown 
515bbd56ea5SKarl Rupp       if ((im < (PetscScalar) is - 1.) || (im > (PetscScalar) ie) || (jm < (PetscScalar)  js - 1.) || (jm > (PetscScalar) je)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "Nonlocal point: (%g,%g)", im, jm);
516af33a6ddSJed Brown     }
517af33a6ddSJed Brown 
518bbd56ea5SKarl Rupp     if (c->fieldInterpLocal) c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
519bbd56ea5SKarl Rupp     else c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
520bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
521af33a6ddSJed Brown   }
522af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
523af33a6ddSJed Brown 
524af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
525af33a6ddSJed Brown   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
526af33a6ddSJed Brown   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
527af33a6ddSJed Brown   if (c->fieldInterpLocal) {
528af33a6ddSJed Brown     ierr = DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
529af33a6ddSJed Brown     ierr = DMRestoreLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
530af33a6ddSJed Brown   }
531af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
532af33a6ddSJed Brown 
533af33a6ddSJed Brown   /* Return field of characteristics at t_n-1 */
534af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
535af33a6ddSJed Brown   ierr = DMDAGetInfo(c->fieldDA,0,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
536af33a6ddSJed Brown   ierr = DMDAVecGetArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
537af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
538af33a6ddSJed Brown     Qi = c->queue[n];
539bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp];
540af33a6ddSJed Brown   }
541af33a6ddSJed Brown   ierr = DMDAVecRestoreArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
542af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
543af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
544af33a6ddSJed Brown 
545af33a6ddSJed Brown   /* Cleanup */
546af33a6ddSJed Brown   ierr = PetscFree(interpIndices);CHKERRQ(ierr);
547af33a6ddSJed Brown   ierr = PetscFree(velocityValues);CHKERRQ(ierr);
548af33a6ddSJed Brown   ierr = PetscFree(velocityValuesOld);CHKERRQ(ierr);
549af33a6ddSJed Brown   ierr = PetscFree(fieldValues);CHKERRQ(ierr);
550af33a6ddSJed Brown   PetscFunctionReturn(0);
551af33a6ddSJed Brown }
552af33a6ddSJed Brown 
553af33a6ddSJed Brown #undef __FUNCT__
554af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetNeighbors"
555af33a6ddSJed Brown PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
556af33a6ddSJed Brown {
557af33a6ddSJed Brown   PetscErrorCode ierr;
558af33a6ddSJed Brown 
559af33a6ddSJed Brown   PetscFunctionBegin;
560af33a6ddSJed Brown   c->numNeighbors = numNeighbors;
5616eaa7dc8SBarry Smith   ierr = PetscFree(c->neighbors);CHKERRQ(ierr);
562af33a6ddSJed Brown   ierr = PetscMalloc(numNeighbors * sizeof(PetscMPIInt), &c->neighbors);CHKERRQ(ierr);
563af33a6ddSJed Brown   ierr = PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));CHKERRQ(ierr);
564af33a6ddSJed Brown   PetscFunctionReturn(0);
565af33a6ddSJed Brown }
566af33a6ddSJed Brown 
567af33a6ddSJed Brown #undef __FUNCT__
568af33a6ddSJed Brown #define __FUNCT__ "CharacteristicAddPoint"
569af33a6ddSJed Brown PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
570af33a6ddSJed Brown {
571af33a6ddSJed Brown   PetscFunctionBegin;
572f23aa3ddSBarry Smith   if (c->queueSize >= c->queueMax) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax);
573af33a6ddSJed Brown   c->queue[c->queueSize++] = *point;
574af33a6ddSJed Brown   PetscFunctionReturn(0);
575af33a6ddSJed Brown }
576af33a6ddSJed Brown 
577af33a6ddSJed Brown #undef __FUNCT__
578af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSendCoordinatesBegin"
579af33a6ddSJed Brown int CharacteristicSendCoordinatesBegin(Characteristic c)
580af33a6ddSJed Brown {
581af33a6ddSJed Brown   PetscMPIInt    rank, tag = 121;
582af33a6ddSJed Brown   PetscInt       i, n;
583af33a6ddSJed Brown   PetscErrorCode ierr;
584af33a6ddSJed Brown 
585af33a6ddSJed Brown   PetscFunctionBegin;
586ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
5870b8f38c8SBarry Smith   ierr = CharacteristicHeapSort(c, c->queue, c->queueSize);CHKERRQ(ierr);
588af33a6ddSJed Brown   ierr = PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));CHKERRQ(ierr);
589bbd56ea5SKarl Rupp   for (i = 0;  i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
590af33a6ddSJed Brown   c->fillCount[0] = 0;
591af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
592ce94432eSBarry Smith     ierr = MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr);
593af33a6ddSJed Brown   }
594af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
595ce94432eSBarry Smith     ierr = MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
596af33a6ddSJed Brown   }
597af33a6ddSJed Brown   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
598af33a6ddSJed Brown   /* Initialize the remote queue */
599af33a6ddSJed Brown   c->queueLocalMax  = c->localOffsets[0]  = 0;
600af33a6ddSJed Brown   c->queueRemoteMax = c->remoteOffsets[0] = 0;
601af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
602af33a6ddSJed Brown     c->remoteOffsets[n] = c->queueRemoteMax;
603af33a6ddSJed Brown     c->queueRemoteMax  += c->fillCount[n];
604af33a6ddSJed Brown     c->localOffsets[n]  = c->queueLocalMax;
605af33a6ddSJed Brown     c->queueLocalMax   += c->needCount[n];
606af33a6ddSJed Brown   }
607af33a6ddSJed Brown   /* HACK BEGIN */
608bbd56ea5SKarl Rupp   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
609af33a6ddSJed Brown   c->needCount[0] = 0;
610af33a6ddSJed Brown   /* HACK END */
611af33a6ddSJed Brown   if (c->queueRemoteMax) {
612af33a6ddSJed Brown     ierr = PetscMalloc(sizeof(CharacteristicPointDA2D) * c->queueRemoteMax, &c->queueRemote);CHKERRQ(ierr);
6130298fd71SBarry Smith   } else c->queueRemote = NULL;
614af33a6ddSJed Brown   c->queueRemoteSize = c->queueRemoteMax;
615af33a6ddSJed Brown 
616af33a6ddSJed Brown   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
617af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
6180298fd71SBarry Smith     ierr = PetscInfo2(NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);CHKERRQ(ierr);
619ce94432eSBarry Smith     ierr = MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr);
620af33a6ddSJed Brown   }
621af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
6220298fd71SBarry Smith     ierr = PetscInfo2(NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);CHKERRQ(ierr);
623ce94432eSBarry Smith     ierr = MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
624af33a6ddSJed Brown   }
625af33a6ddSJed Brown   PetscFunctionReturn(0);
626af33a6ddSJed Brown }
627af33a6ddSJed Brown 
628af33a6ddSJed Brown #undef __FUNCT__
629af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSendCoordinatesEnd"
630af33a6ddSJed Brown PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
631af33a6ddSJed Brown {
632af33a6ddSJed Brown #if 0
633af33a6ddSJed Brown   PetscMPIInt rank;
634af33a6ddSJed Brown   PetscInt    n;
635af33a6ddSJed Brown #endif
636af33a6ddSJed Brown   PetscErrorCode ierr;
637af33a6ddSJed Brown 
638af33a6ddSJed Brown   PetscFunctionBegin;
639af33a6ddSJed Brown   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
640af33a6ddSJed Brown #if 0
641ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
642af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
643f23aa3ddSBarry Smith     if (c->neighbors[c->queueRemote[n].proc] == rank) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "This is messed up, n = %d proc = %d", n, c->queueRemote[n].proc);
644af33a6ddSJed Brown   }
645af33a6ddSJed Brown #endif
646af33a6ddSJed Brown   PetscFunctionReturn(0);
647af33a6ddSJed Brown }
648af33a6ddSJed Brown 
649af33a6ddSJed Brown #undef __FUNCT__
650af33a6ddSJed Brown #define __FUNCT__ "CharacteristicGetValuesBegin"
651af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
652af33a6ddSJed Brown {
653af33a6ddSJed Brown   PetscMPIInt    tag = 121;
654af33a6ddSJed Brown   PetscInt       n;
655af33a6ddSJed Brown   PetscErrorCode ierr;
656af33a6ddSJed Brown 
657af33a6ddSJed Brown   PetscFunctionBegin;
658af33a6ddSJed Brown   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
659af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
660ce94432eSBarry Smith     ierr = MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr);
661af33a6ddSJed Brown   }
662af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
663ce94432eSBarry Smith     ierr = MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
664af33a6ddSJed Brown   }
665af33a6ddSJed Brown   PetscFunctionReturn(0);
666af33a6ddSJed Brown }
667af33a6ddSJed Brown 
668af33a6ddSJed Brown #undef __FUNCT__
669af33a6ddSJed Brown #define __FUNCT__ "CharacteristicGetValuesEnd"
670af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
671af33a6ddSJed Brown {
672af33a6ddSJed Brown   PetscErrorCode ierr;
673af33a6ddSJed Brown 
674af33a6ddSJed Brown   PetscFunctionBegin;
675af33a6ddSJed Brown   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
676af33a6ddSJed Brown   /* Free queue of requests from other procs */
677af33a6ddSJed Brown   ierr = PetscFree(c->queueRemote);CHKERRQ(ierr);
678af33a6ddSJed Brown   PetscFunctionReturn(0);
679af33a6ddSJed Brown }
680af33a6ddSJed Brown 
681af33a6ddSJed Brown /*---------------------------------------------------------------------*/
682af33a6ddSJed Brown #undef __FUNCT__
6830b8f38c8SBarry Smith #define __FUNCT__ "CharacteristicHeapSort"
684af33a6ddSJed Brown /*
685af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
686af33a6ddSJed Brown */
6870b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
688af33a6ddSJed Brown /*---------------------------------------------------------------------*/
689af33a6ddSJed Brown {
6900b8f38c8SBarry Smith   PetscErrorCode          ierr;
691af33a6ddSJed Brown   CharacteristicPointDA2D temp;
6920b8f38c8SBarry Smith   PetscInt                n;
693af33a6ddSJed Brown 
6940b8f38c8SBarry Smith   PetscFunctionBegin;
695af33a6ddSJed Brown   if (0) { /* Check the order of the queue before sorting */
6960298fd71SBarry Smith     PetscInfo(NULL, "Before Heap sort\n");
697af33a6ddSJed Brown     for (n=0; n<size; n++) {
6980298fd71SBarry Smith       ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
699af33a6ddSJed Brown     }
700af33a6ddSJed Brown   }
701af33a6ddSJed Brown 
702af33a6ddSJed Brown   /* SORTING PHASE */
7030b8f38c8SBarry Smith   for (n = (size / 2)-1; n >= 0; n--) {
7040b8f38c8SBarry Smith     ierr = CharacteristicSiftDown(c, queue, n, size-1);CHKERRQ(ierr); /* Rich had size-1 here, Matt had size*/
7050b8f38c8SBarry Smith   }
706af33a6ddSJed Brown   for (n = size-1; n >= 1; n--) {
707af33a6ddSJed Brown     temp     = queue[0];
708af33a6ddSJed Brown     queue[0] = queue[n];
709af33a6ddSJed Brown     queue[n] = temp;
7100b8f38c8SBarry Smith     ierr     = CharacteristicSiftDown(c, queue, 0, n-1);CHKERRQ(ierr);
711af33a6ddSJed Brown   }
712af33a6ddSJed Brown   if (0) { /* Check the order of the queue after sorting */
7130298fd71SBarry Smith     ierr = PetscInfo(NULL, "Avter  Heap sort\n");CHKERRQ(ierr);
714af33a6ddSJed Brown     for (n=0; n<size; n++) {
7150298fd71SBarry Smith       ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
716af33a6ddSJed Brown     }
717af33a6ddSJed Brown   }
7180b8f38c8SBarry Smith   PetscFunctionReturn(0);
719af33a6ddSJed Brown }
720af33a6ddSJed Brown 
721af33a6ddSJed Brown /*---------------------------------------------------------------------*/
722af33a6ddSJed Brown #undef __FUNCT__
7230b8f38c8SBarry Smith #define __FUNCT__ "CharacteristicSiftDown"
724af33a6ddSJed Brown /*
725af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
726af33a6ddSJed Brown */
7270b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
728af33a6ddSJed Brown /*---------------------------------------------------------------------*/
729af33a6ddSJed Brown {
730af33a6ddSJed Brown   PetscBool               done = PETSC_FALSE;
731af33a6ddSJed Brown   PetscInt                maxChild;
732af33a6ddSJed Brown   CharacteristicPointDA2D temp;
733af33a6ddSJed Brown 
734af33a6ddSJed Brown   PetscFunctionBegin;
735af33a6ddSJed Brown   while ((root*2 <= bottom) && (!done)) {
736af33a6ddSJed Brown     if (root*2 == bottom) maxChild = root * 2;
737af33a6ddSJed Brown     else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2;
738af33a6ddSJed Brown     else maxChild = root * 2 + 1;
739af33a6ddSJed Brown 
740af33a6ddSJed Brown     if (queue[root].proc < queue[maxChild].proc) {
741af33a6ddSJed Brown       temp            = queue[root];
742af33a6ddSJed Brown       queue[root]     = queue[maxChild];
743af33a6ddSJed Brown       queue[maxChild] = temp;
744af33a6ddSJed Brown       root            = maxChild;
745db4deed7SKarl Rupp     } else done = PETSC_TRUE;
746af33a6ddSJed Brown   }
747af33a6ddSJed Brown   PetscFunctionReturn(0);
748af33a6ddSJed Brown }
749af33a6ddSJed Brown 
750af33a6ddSJed Brown #undef __FUNCT__
751af33a6ddSJed Brown #define __FUNCT__ "DMDAGetNeighborsRank"
752af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
753af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
754af33a6ddSJed Brown {
755af33a6ddSJed Brown   DMDABoundaryType bx, by;
756af33a6ddSJed Brown   PetscBool        IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
757af33a6ddSJed Brown   MPI_Comm         comm;
758af33a6ddSJed Brown   PetscMPIInt      rank;
759af33a6ddSJed Brown   PetscInt         **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ;
760af33a6ddSJed Brown   PetscErrorCode   ierr;
761af33a6ddSJed Brown 
762af33a6ddSJed Brown   PetscFunctionBegin;
763af33a6ddSJed Brown   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
764af33a6ddSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
76522d28d08SBarry Smith   ierr = DMDAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &bx, &by,0, 0);CHKERRQ(ierr);
766af33a6ddSJed Brown 
767bbd56ea5SKarl Rupp   if (bx == DMDA_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
768bbd56ea5SKarl Rupp   if (by == DMDA_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
769af33a6ddSJed Brown 
770af33a6ddSJed Brown   neighbors[0] = rank;
771af33a6ddSJed Brown   rank = 0;
772af33a6ddSJed Brown   ierr = PetscMalloc(sizeof(PetscInt*)*PJ,&procs);CHKERRQ(ierr);
773af33a6ddSJed Brown   for (pj=0; pj<PJ; pj++) {
774af33a6ddSJed Brown     ierr = PetscMalloc(sizeof(PetscInt)*PI,&(procs[pj]));CHKERRQ(ierr);
775af33a6ddSJed Brown     for (pi=0; pi<PI; pi++) {
776af33a6ddSJed Brown       procs[pj][pi] = rank;
777af33a6ddSJed Brown       rank++;
778af33a6ddSJed Brown     }
779af33a6ddSJed Brown   }
780af33a6ddSJed Brown 
781af33a6ddSJed Brown   pi  = neighbors[0] % PI;
782af33a6ddSJed Brown   pj  = neighbors[0] / PI;
783af33a6ddSJed Brown   pim = pi-1;  if (pim<0) pim=PI-1;
784af33a6ddSJed Brown   pip = (pi+1)%PI;
785af33a6ddSJed Brown   pjm = pj-1;  if (pjm<0) pjm=PJ-1;
786af33a6ddSJed Brown   pjp = (pj+1)%PJ;
787af33a6ddSJed Brown 
788af33a6ddSJed Brown   neighbors[1] = procs[pj] [pim];
789af33a6ddSJed Brown   neighbors[2] = procs[pjp][pim];
790af33a6ddSJed Brown   neighbors[3] = procs[pjp][pi];
791af33a6ddSJed Brown   neighbors[4] = procs[pjp][pip];
792af33a6ddSJed Brown   neighbors[5] = procs[pj] [pip];
793af33a6ddSJed Brown   neighbors[6] = procs[pjm][pip];
794af33a6ddSJed Brown   neighbors[7] = procs[pjm][pi];
795af33a6ddSJed Brown   neighbors[8] = procs[pjm][pim];
796af33a6ddSJed Brown 
797af33a6ddSJed Brown   if (!IPeriodic) {
798af33a6ddSJed Brown     if (pi==0)    neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0];
799af33a6ddSJed Brown     if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0];
800af33a6ddSJed Brown   }
801af33a6ddSJed Brown 
802af33a6ddSJed Brown   if (!JPeriodic) {
803af33a6ddSJed Brown     if (pj==0)    neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0];
804af33a6ddSJed Brown     if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0];
805af33a6ddSJed Brown   }
806af33a6ddSJed Brown 
807af33a6ddSJed Brown   for (pj = 0; pj < PJ; pj++) {
808af33a6ddSJed Brown     ierr = PetscFree(procs[pj]);CHKERRQ(ierr);
809af33a6ddSJed Brown   }
810af33a6ddSJed Brown   ierr = PetscFree(procs);CHKERRQ(ierr);
811af33a6ddSJed Brown   PetscFunctionReturn(0);
812af33a6ddSJed Brown }
813af33a6ddSJed Brown 
814af33a6ddSJed Brown #undef __FUNCT__
815af33a6ddSJed Brown #define __FUNCT__ "DMDAGetNeighborRelative"
816af33a6ddSJed Brown /*
817af33a6ddSJed Brown   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
818af33a6ddSJed Brown     2 | 3 | 4
819af33a6ddSJed Brown     __|___|__
820af33a6ddSJed Brown     1 | 0 | 5
821af33a6ddSJed Brown     __|___|__
822af33a6ddSJed Brown     8 | 7 | 6
823af33a6ddSJed Brown       |   |
824af33a6ddSJed Brown */
825af33a6ddSJed Brown PetscInt DMDAGetNeighborRelative(DM da, PassiveReal ir, PassiveReal jr)
826af33a6ddSJed Brown {
827af33a6ddSJed Brown   DMDALocalInfo  info;
828af33a6ddSJed Brown   PassiveReal    is,ie,js,je;
829af33a6ddSJed Brown   PetscErrorCode ierr;
830af33a6ddSJed Brown 
831af33a6ddSJed Brown   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
832af33a6ddSJed Brown   is   = (PassiveReal) info.xs - 0.5; ie = (PassiveReal) info.xs + info.xm - 0.5;
833af33a6ddSJed Brown   js   = (PassiveReal) info.ys - 0.5; je = (PassiveReal) info.ys + info.ym - 0.5;
834af33a6ddSJed Brown 
835af33a6ddSJed Brown   if (ir >= is && ir <= ie) { /* center column */
836bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 0;
837bbd56ea5SKarl Rupp     else if (jr < js)         return 7;
838bbd56ea5SKarl Rupp     else                      return 3;
839af33a6ddSJed Brown   } else if (ir < is) {     /* left column */
840bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 1;
841bbd56ea5SKarl Rupp     else if (jr < js)         return 8;
842bbd56ea5SKarl Rupp     else                      return 2;
843af33a6ddSJed Brown   } else {                  /* right column */
844bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 5;
845bbd56ea5SKarl Rupp     else if (jr < js)         return 6;
846bbd56ea5SKarl Rupp     else                      return 4;
847af33a6ddSJed Brown   }
848af33a6ddSJed Brown }
849