xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision 61ba467669c3da7ee29950d5012260364c1524d5)
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 #undef __FUNCT__
24af33a6ddSJed Brown #define __FUNCT__ "CharacteristicView"
25af33a6ddSJed Brown PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
26af33a6ddSJed Brown {
27af33a6ddSJed Brown   PetscBool      iascii;
28af33a6ddSJed Brown   PetscErrorCode ierr;
29af33a6ddSJed Brown 
30af33a6ddSJed Brown   PetscFunctionBegin;
31af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
32af33a6ddSJed Brown   if (!viewer) {
33ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr);
34af33a6ddSJed Brown   }
35af33a6ddSJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
36af33a6ddSJed Brown   PetscCheckSameComm(c, 1, viewer, 2);
37af33a6ddSJed Brown 
38251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
39af33a6ddSJed Brown   if (!iascii) {
40af33a6ddSJed Brown     if (c->ops->view) {
41af33a6ddSJed Brown       ierr = (*c->ops->view)(c, viewer);CHKERRQ(ierr);
42af33a6ddSJed Brown     }
43af33a6ddSJed Brown   }
44af33a6ddSJed Brown   PetscFunctionReturn(0);
45af33a6ddSJed Brown }
46af33a6ddSJed Brown 
47af33a6ddSJed Brown #undef __FUNCT__
48af33a6ddSJed Brown #define __FUNCT__ "CharacteristicDestroy"
496bf464f9SBarry Smith PetscErrorCode CharacteristicDestroy(Characteristic *c)
50af33a6ddSJed Brown {
51af33a6ddSJed Brown   PetscErrorCode ierr;
52af33a6ddSJed Brown 
53af33a6ddSJed Brown   PetscFunctionBegin;
546bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
556bf464f9SBarry Smith   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
566bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0);
57af33a6ddSJed Brown 
586bf464f9SBarry Smith   if ((*c)->ops->destroy) {
596bf464f9SBarry Smith     ierr = (*(*c)->ops->destroy)((*c));CHKERRQ(ierr);
60af33a6ddSJed Brown   }
616bf464f9SBarry Smith   ierr = MPI_Type_free(&(*c)->itemType);CHKERRQ(ierr);
626bf464f9SBarry Smith   ierr = PetscFree((*c)->queue);CHKERRQ(ierr);
636bf464f9SBarry Smith   ierr = PetscFree((*c)->queueLocal);CHKERRQ(ierr);
646bf464f9SBarry Smith   ierr = PetscFree((*c)->queueRemote);CHKERRQ(ierr);
656bf464f9SBarry Smith   ierr = PetscFree((*c)->neighbors);CHKERRQ(ierr);
666bf464f9SBarry Smith   ierr = PetscFree((*c)->needCount);CHKERRQ(ierr);
676bf464f9SBarry Smith   ierr = PetscFree((*c)->localOffsets);CHKERRQ(ierr);
686bf464f9SBarry Smith   ierr = PetscFree((*c)->fillCount);CHKERRQ(ierr);
696bf464f9SBarry Smith   ierr = PetscFree((*c)->remoteOffsets);CHKERRQ(ierr);
706bf464f9SBarry Smith   ierr = PetscFree((*c)->request);CHKERRQ(ierr);
716bf464f9SBarry Smith   ierr = PetscFree((*c)->status);CHKERRQ(ierr);
72af33a6ddSJed Brown   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
73af33a6ddSJed Brown   PetscFunctionReturn(0);
74af33a6ddSJed Brown }
75af33a6ddSJed Brown 
76af33a6ddSJed Brown #undef __FUNCT__
77af33a6ddSJed Brown #define __FUNCT__ "CharacteristicCreate"
78af33a6ddSJed Brown PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
79af33a6ddSJed Brown {
80af33a6ddSJed Brown   Characteristic newC;
81af33a6ddSJed Brown   PetscErrorCode ierr;
82af33a6ddSJed Brown 
83af33a6ddSJed Brown   PetscFunctionBegin;
84af33a6ddSJed Brown   PetscValidPointer(c, 2);
850298fd71SBarry Smith   *c = NULL;
86607a6623SBarry Smith   ierr = CharacteristicInitializePackage();CHKERRQ(ierr);
87af33a6ddSJed Brown 
8873107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "SemiLagrange", comm, CharacteristicDestroy, CharacteristicView);CHKERRQ(ierr);
89af33a6ddSJed Brown   *c   = newC;
90af33a6ddSJed Brown 
91af33a6ddSJed Brown   newC->structured          = PETSC_TRUE;
92af33a6ddSJed Brown   newC->numIds              = 0;
930298fd71SBarry Smith   newC->velocityDA          = NULL;
940298fd71SBarry Smith   newC->velocity            = NULL;
950298fd71SBarry Smith   newC->velocityOld         = NULL;
96af33a6ddSJed Brown   newC->numVelocityComp     = 0;
970298fd71SBarry Smith   newC->velocityComp        = NULL;
980298fd71SBarry Smith   newC->velocityInterp      = NULL;
990298fd71SBarry Smith   newC->velocityInterpLocal = NULL;
1000298fd71SBarry Smith   newC->velocityCtx         = NULL;
1010298fd71SBarry Smith   newC->fieldDA             = NULL;
1020298fd71SBarry Smith   newC->field               = NULL;
103af33a6ddSJed Brown   newC->numFieldComp        = 0;
1040298fd71SBarry Smith   newC->fieldComp           = NULL;
1050298fd71SBarry Smith   newC->fieldInterp         = NULL;
1060298fd71SBarry Smith   newC->fieldInterpLocal    = NULL;
1070298fd71SBarry Smith   newC->fieldCtx            = NULL;
1080298fd71SBarry Smith   newC->itemType            = 0;
1090298fd71SBarry Smith   newC->queue               = NULL;
110af33a6ddSJed Brown   newC->queueSize           = 0;
111af33a6ddSJed Brown   newC->queueMax            = 0;
1120298fd71SBarry Smith   newC->queueLocal          = NULL;
113af33a6ddSJed Brown   newC->queueLocalSize      = 0;
114af33a6ddSJed Brown   newC->queueLocalMax       = 0;
1150298fd71SBarry Smith   newC->queueRemote         = NULL;
116af33a6ddSJed Brown   newC->queueRemoteSize     = 0;
117af33a6ddSJed Brown   newC->queueRemoteMax      = 0;
118af33a6ddSJed Brown   newC->numNeighbors        = 0;
1190298fd71SBarry Smith   newC->neighbors           = NULL;
1200298fd71SBarry Smith   newC->needCount           = NULL;
1210298fd71SBarry Smith   newC->localOffsets        = NULL;
1220298fd71SBarry Smith   newC->fillCount           = NULL;
1230298fd71SBarry Smith   newC->remoteOffsets       = NULL;
1240298fd71SBarry Smith   newC->request             = NULL;
1250298fd71SBarry Smith   newC->status              = NULL;
126af33a6ddSJed Brown   PetscFunctionReturn(0);
127af33a6ddSJed Brown }
128af33a6ddSJed Brown 
129af33a6ddSJed Brown #undef __FUNCT__
130af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetType"
131af33a6ddSJed Brown /*@C
132af33a6ddSJed Brown    CharacteristicSetType - Builds Characteristic for a particular solver.
133af33a6ddSJed Brown 
134af33a6ddSJed Brown    Logically Collective on Characteristic
135af33a6ddSJed Brown 
136af33a6ddSJed Brown    Input Parameters:
137af33a6ddSJed Brown +  c    - the method of characteristics context
138af33a6ddSJed Brown -  type - a known method
139af33a6ddSJed Brown 
140af33a6ddSJed Brown    Options Database Key:
141af33a6ddSJed Brown .  -characteristic_type <method> - Sets the method; use -help for a list
142af33a6ddSJed Brown     of available methods
143af33a6ddSJed Brown 
144af33a6ddSJed Brown    Notes:
14530a76a96SBarry Smith    See "include/petsccharacteristic.h" for available methods
146af33a6ddSJed Brown 
147af33a6ddSJed Brown   Normally, it is best to use the CharacteristicSetFromOptions() command and
148af33a6ddSJed Brown   then set the Characteristic type from the options database rather than by using
149af33a6ddSJed Brown   this routine.  Using the options database provides the user with
150af33a6ddSJed Brown   maximum flexibility in evaluating the many different Krylov methods.
151af33a6ddSJed Brown   The CharacteristicSetType() routine is provided for those situations where it
152af33a6ddSJed Brown   is necessary to set the iterative solver independently of the command
153af33a6ddSJed Brown   line or options database.  This might be the case, for example, when
154af33a6ddSJed Brown   the choice of iterative solver changes during the execution of the
155af33a6ddSJed Brown   program, and the user's application is taking responsibility for
156af33a6ddSJed Brown   choosing the appropriate method.  In other words, this routine is
157af33a6ddSJed Brown   not for beginners.
158af33a6ddSJed Brown 
159af33a6ddSJed Brown   Level: intermediate
160af33a6ddSJed Brown 
161af33a6ddSJed Brown .keywords: Characteristic, set, method
162af33a6ddSJed Brown 
163af33a6ddSJed Brown .seealso: CharacteristicType
164af33a6ddSJed Brown 
165af33a6ddSJed Brown @*/
16619fd82e9SBarry Smith PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
167af33a6ddSJed Brown {
168af33a6ddSJed Brown   PetscErrorCode ierr, (*r)(Characteristic);
169af33a6ddSJed Brown   PetscBool      match;
170af33a6ddSJed Brown 
171af33a6ddSJed Brown   PetscFunctionBegin;
172af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
173af33a6ddSJed Brown   PetscValidCharPointer(type, 2);
174af33a6ddSJed Brown 
175251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) c, type, &match);CHKERRQ(ierr);
176af33a6ddSJed Brown   if (match) PetscFunctionReturn(0);
177af33a6ddSJed Brown 
178af33a6ddSJed Brown   if (c->data) {
179af33a6ddSJed Brown     /* destroy the old private Characteristic context */
180af33a6ddSJed Brown     ierr = (*c->ops->destroy)(c);CHKERRQ(ierr);
1810298fd71SBarry Smith     c->ops->destroy = NULL;
182af33a6ddSJed Brown     c->data         = 0;
183af33a6ddSJed Brown   }
184af33a6ddSJed Brown 
1851c9cd337SJed Brown   ierr =  PetscFunctionListFind(CharacteristicList,type,&r);CHKERRQ(ierr);
186af33a6ddSJed Brown   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
187af33a6ddSJed Brown   c->setupcalled = 0;
188af33a6ddSJed Brown   ierr = (*r)(c);CHKERRQ(ierr);
189af33a6ddSJed Brown   ierr = PetscObjectChangeTypeName((PetscObject) c, type);CHKERRQ(ierr);
190af33a6ddSJed Brown   PetscFunctionReturn(0);
191af33a6ddSJed Brown }
192af33a6ddSJed Brown 
193af33a6ddSJed Brown #undef __FUNCT__
194af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetUp"
195af33a6ddSJed Brown /*@
196af33a6ddSJed Brown    CharacteristicSetUp - Sets up the internal data structures for the
197af33a6ddSJed Brown    later use of an iterative solver.
198af33a6ddSJed Brown 
199af33a6ddSJed Brown    Collective on Characteristic
200af33a6ddSJed Brown 
201af33a6ddSJed Brown    Input Parameter:
202af33a6ddSJed Brown .  ksp   - iterative context obtained from CharacteristicCreate()
203af33a6ddSJed Brown 
204af33a6ddSJed Brown    Level: developer
205af33a6ddSJed Brown 
206af33a6ddSJed Brown .keywords: Characteristic, setup
207af33a6ddSJed Brown 
208af33a6ddSJed Brown .seealso: CharacteristicCreate(), CharacteristicSolve(), CharacteristicDestroy()
209af33a6ddSJed Brown @*/
210af33a6ddSJed Brown PetscErrorCode CharacteristicSetUp(Characteristic c)
211af33a6ddSJed Brown {
212af33a6ddSJed Brown   PetscErrorCode ierr;
213af33a6ddSJed Brown 
214af33a6ddSJed Brown   PetscFunctionBegin;
215af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
216af33a6ddSJed Brown 
217af33a6ddSJed Brown   if (!((PetscObject)c)->type_name) {
218af33a6ddSJed Brown     ierr = CharacteristicSetType(c, CHARACTERISTICDA);CHKERRQ(ierr);
219af33a6ddSJed Brown   }
220af33a6ddSJed Brown 
221af33a6ddSJed Brown   if (c->setupcalled == 2) PetscFunctionReturn(0);
222af33a6ddSJed Brown 
223af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr);
224af33a6ddSJed Brown   if (!c->setupcalled) {
225af33a6ddSJed Brown     ierr = (*c->ops->setup)(c);CHKERRQ(ierr);
226af33a6ddSJed Brown   }
227af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr);
228af33a6ddSJed Brown   c->setupcalled = 2;
229af33a6ddSJed Brown   PetscFunctionReturn(0);
230af33a6ddSJed Brown }
231af33a6ddSJed Brown 
232af33a6ddSJed Brown #undef __FUNCT__
233af33a6ddSJed Brown #define __FUNCT__ "CharacteristicRegister"
234af33a6ddSJed Brown /*@C
2351c84c290SBarry Smith   CharacteristicRegister -  Adds a solver to the method of characteristics package.
2361c84c290SBarry Smith 
2371c84c290SBarry Smith    Not Collective
2381c84c290SBarry Smith 
2391c84c290SBarry Smith    Input Parameters:
2401c84c290SBarry Smith +  name_solver - name of a new user-defined solver
2411c84c290SBarry Smith -  routine_create - routine to create method context
2421c84c290SBarry Smith 
2431c84c290SBarry Smith   Sample usage:
2441c84c290SBarry Smith .vb
245bdf89e91SBarry Smith     CharacteristicRegister("my_char", MyCharCreate);
2461c84c290SBarry Smith .ve
2471c84c290SBarry Smith 
2481c84c290SBarry Smith   Then, your Characteristic type can be chosen with the procedural interface via
2491c84c290SBarry Smith .vb
2501c84c290SBarry Smith     CharacteristicCreate(MPI_Comm, Characteristic* &char);
2511c84c290SBarry Smith     CharacteristicSetType(char,"my_char");
2521c84c290SBarry Smith .ve
2531c84c290SBarry Smith    or at runtime via the option
2541c84c290SBarry Smith .vb
2551c84c290SBarry Smith     -characteristic_type my_char
2561c84c290SBarry Smith .ve
2571c84c290SBarry Smith 
2581c84c290SBarry Smith    Notes:
2591c84c290SBarry Smith    CharacteristicRegister() may be called multiple times to add several user-defined solvers.
2601c84c290SBarry Smith 
2611c84c290SBarry Smith .keywords: Characteristic, register
2621c84c290SBarry Smith 
2631c84c290SBarry Smith .seealso: CharacteristicRegisterAll(), CharacteristicRegisterDestroy()
264af33a6ddSJed Brown 
265af33a6ddSJed Brown   Level: advanced
266af33a6ddSJed Brown @*/
267bdf89e91SBarry Smith PetscErrorCode CharacteristicRegister(const char sname[],PetscErrorCode (*function)(Characteristic))
268af33a6ddSJed Brown {
269af33a6ddSJed Brown   PetscErrorCode ierr;
270af33a6ddSJed Brown 
271af33a6ddSJed Brown   PetscFunctionBegin;
272a240a19fSJed Brown   ierr = PetscFunctionListAdd(&CharacteristicList,sname,function);CHKERRQ(ierr);
273af33a6ddSJed Brown   PetscFunctionReturn(0);
274af33a6ddSJed Brown }
275af33a6ddSJed Brown 
276af33a6ddSJed Brown #undef __FUNCT__
277af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetVelocityInterpolation"
278af33a6ddSJed 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)
279af33a6ddSJed Brown {
280af33a6ddSJed Brown   PetscFunctionBegin;
281af33a6ddSJed Brown   c->velocityDA      = da;
282af33a6ddSJed Brown   c->velocity        = v;
283af33a6ddSJed Brown   c->velocityOld     = vOld;
284af33a6ddSJed Brown   c->numVelocityComp = numComponents;
285af33a6ddSJed Brown   c->velocityComp    = components;
286af33a6ddSJed Brown   c->velocityInterp  = interp;
287af33a6ddSJed Brown   c->velocityCtx     = ctx;
288af33a6ddSJed Brown   PetscFunctionReturn(0);
289af33a6ddSJed Brown }
290af33a6ddSJed Brown 
291af33a6ddSJed Brown #undef __FUNCT__
292af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetVelocityInterpolationLocal"
293af33a6ddSJed 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)
294af33a6ddSJed Brown {
295af33a6ddSJed Brown   PetscFunctionBegin;
296af33a6ddSJed Brown   c->velocityDA          = da;
297af33a6ddSJed Brown   c->velocity            = v;
298af33a6ddSJed Brown   c->velocityOld         = vOld;
299af33a6ddSJed Brown   c->numVelocityComp     = numComponents;
300af33a6ddSJed Brown   c->velocityComp        = components;
301af33a6ddSJed Brown   c->velocityInterpLocal = interp;
302af33a6ddSJed Brown   c->velocityCtx         = ctx;
303af33a6ddSJed Brown   PetscFunctionReturn(0);
304af33a6ddSJed Brown }
305af33a6ddSJed Brown 
306af33a6ddSJed Brown #undef __FUNCT__
307af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetFieldInterpolation"
308af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
309af33a6ddSJed Brown {
310af33a6ddSJed Brown   PetscFunctionBegin;
311af33a6ddSJed Brown #if 0
312af33a6ddSJed 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.");
313af33a6ddSJed Brown #endif
314af33a6ddSJed Brown   c->fieldDA      = da;
315af33a6ddSJed Brown   c->field        = v;
316af33a6ddSJed Brown   c->numFieldComp = numComponents;
317af33a6ddSJed Brown   c->fieldComp    = components;
318af33a6ddSJed Brown   c->fieldInterp  = interp;
319af33a6ddSJed Brown   c->fieldCtx     = ctx;
320af33a6ddSJed Brown   PetscFunctionReturn(0);
321af33a6ddSJed Brown }
322af33a6ddSJed Brown 
323af33a6ddSJed Brown #undef __FUNCT__
324af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetFieldInterpolationLocal"
325af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal[], PetscInt, PetscInt[], PetscScalar [], void*), void *ctx)
326af33a6ddSJed Brown {
327af33a6ddSJed Brown   PetscFunctionBegin;
328af33a6ddSJed Brown #if 0
329af33a6ddSJed 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.");
330af33a6ddSJed Brown #endif
331af33a6ddSJed Brown   c->fieldDA          = da;
332af33a6ddSJed Brown   c->field            = v;
333af33a6ddSJed Brown   c->numFieldComp     = numComponents;
334af33a6ddSJed Brown   c->fieldComp        = components;
335af33a6ddSJed Brown   c->fieldInterpLocal = interp;
336af33a6ddSJed Brown   c->fieldCtx         = ctx;
337af33a6ddSJed Brown   PetscFunctionReturn(0);
338af33a6ddSJed Brown }
339af33a6ddSJed Brown 
340af33a6ddSJed Brown #undef __FUNCT__
341af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSolve"
342af33a6ddSJed Brown PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
343af33a6ddSJed Brown {
344af33a6ddSJed Brown   CharacteristicPointDA2D Qi;
345af33a6ddSJed Brown   DM                      da = c->velocityDA;
346af33a6ddSJed Brown   Vec                     velocityLocal, velocityLocalOld;
347af33a6ddSJed Brown   Vec                     fieldLocal;
348af33a6ddSJed Brown   DMDALocalInfo           info;
349af33a6ddSJed Brown   PetscScalar             **solArray;
350af33a6ddSJed Brown   void                    *velocityArray;
351af33a6ddSJed Brown   void                    *velocityArrayOld;
352af33a6ddSJed Brown   void                    *fieldArray;
3531cc8b206SBarry Smith   PetscScalar             *interpIndices;
3541cc8b206SBarry Smith   PetscScalar             *velocityValues, *velocityValuesOld;
3551cc8b206SBarry Smith   PetscScalar             *fieldValues;
356af33a6ddSJed Brown   PetscMPIInt             rank;
357af33a6ddSJed Brown   PetscInt                dim;
358af33a6ddSJed Brown   PetscMPIInt             neighbors[9];
359af33a6ddSJed Brown   PetscInt                dof;
360af33a6ddSJed Brown   PetscInt                gx, gy;
361af33a6ddSJed Brown   PetscInt                n, is, ie, js, je, comp;
362af33a6ddSJed Brown   PetscErrorCode          ierr;
363af33a6ddSJed Brown 
364af33a6ddSJed Brown   PetscFunctionBegin;
365af33a6ddSJed Brown   c->queueSize = 0;
366ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
367af33a6ddSJed Brown   ierr = DMDAGetNeighborsRank(da, neighbors);CHKERRQ(ierr);
368af33a6ddSJed Brown   ierr = CharacteristicSetNeighbors(c, 9, neighbors);CHKERRQ(ierr);
369af33a6ddSJed Brown   ierr = CharacteristicSetUp(c);CHKERRQ(ierr);
370af33a6ddSJed Brown   /* global and local grid info */
371af33a6ddSJed Brown   ierr = DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);CHKERRQ(ierr);
372af33a6ddSJed Brown   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
373af33a6ddSJed Brown   is   = info.xs;          ie   = info.xs+info.xm;
374af33a6ddSJed Brown   js   = info.ys;          je   = info.ys+info.ym;
375af33a6ddSJed Brown   /* Allocation */
376785e854fSJed Brown   ierr = PetscMalloc1(dim,                &interpIndices);CHKERRQ(ierr);
377785e854fSJed Brown   ierr = PetscMalloc1(c->numVelocityComp, &velocityValues);CHKERRQ(ierr);
378785e854fSJed Brown   ierr = PetscMalloc1(c->numVelocityComp, &velocityValuesOld);CHKERRQ(ierr);
379785e854fSJed Brown   ierr = PetscMalloc1(c->numFieldComp,    &fieldValues);CHKERRQ(ierr);
380af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
381af33a6ddSJed Brown 
382af33a6ddSJed Brown   /* -----------------------------------------------------------------------
383af33a6ddSJed Brown      PART 1, AT t-dt/2
384af33a6ddSJed Brown      -----------------------------------------------------------------------*/
385af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
386af33a6ddSJed Brown   /* GET POSITION AT HALF TIME IN THE PAST */
387af33a6ddSJed Brown   if (c->velocityInterpLocal) {
388af33a6ddSJed Brown     ierr = DMGetLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
389af33a6ddSJed Brown     ierr = DMGetLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
390af33a6ddSJed Brown     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
391af33a6ddSJed Brown     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
392af33a6ddSJed Brown     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
393af33a6ddSJed Brown     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
394af33a6ddSJed Brown     ierr = DMDAVecGetArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
395af33a6ddSJed Brown     ierr = DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
396af33a6ddSJed Brown   }
3970298fd71SBarry Smith   ierr = PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n");CHKERRQ(ierr);
398af33a6ddSJed Brown   for (Qi.j = js; Qi.j < je; Qi.j++) {
399af33a6ddSJed Brown     for (Qi.i = is; Qi.i < ie; Qi.i++) {
400af33a6ddSJed Brown       interpIndices[0] = Qi.i;
401af33a6ddSJed Brown       interpIndices[1] = Qi.j;
402*61ba4676SBarry Smith       if (c->velocityInterpLocal) {ierr = c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr);}
403*61ba4676SBarry Smith       else {ierr = c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr);}
404af33a6ddSJed Brown       Qi.x = Qi.i - velocityValues[0]*dt/2.0;
405af33a6ddSJed Brown       Qi.y = Qi.j - velocityValues[1]*dt/2.0;
406af33a6ddSJed Brown 
407af33a6ddSJed Brown       /* Determine whether the position at t - dt/2 is local */
408af33a6ddSJed Brown       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
409af33a6ddSJed Brown 
410af33a6ddSJed Brown       /* Check for Periodic boundaries and move all periodic points back onto the domain */
411af33a6ddSJed Brown       ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
412af33a6ddSJed Brown       ierr = CharacteristicAddPoint(c, &Qi);CHKERRQ(ierr);
413af33a6ddSJed Brown     }
414af33a6ddSJed Brown   }
415af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
416af33a6ddSJed Brown 
417af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
418af33a6ddSJed Brown   ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr);
419af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
420af33a6ddSJed Brown 
421af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
422af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (local values) */
4230298fd71SBarry Smith   ierr = PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n");CHKERRQ(ierr);
424af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
425af33a6ddSJed Brown     Qi = c->queue[n];
426af33a6ddSJed Brown     if (c->neighbors[Qi.proc] == rank) {
427af33a6ddSJed Brown       interpIndices[0] = Qi.x;
428af33a6ddSJed Brown       interpIndices[1] = Qi.y;
429af33a6ddSJed Brown       if (c->velocityInterpLocal) {
430*61ba4676SBarry Smith         ierr = c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr);
431*61ba4676SBarry Smith         ierr = c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);CHKERRQ(ierr);
432af33a6ddSJed Brown       } else {
433*61ba4676SBarry Smith         ierr = c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr);
434*61ba4676SBarry Smith         ierr = c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
435af33a6ddSJed Brown       }
436af33a6ddSJed Brown       Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
437af33a6ddSJed Brown       Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
438af33a6ddSJed Brown     }
439af33a6ddSJed Brown     c->queue[n] = Qi;
440af33a6ddSJed Brown   }
441af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
442af33a6ddSJed Brown 
443af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
444af33a6ddSJed Brown   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
445af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
446af33a6ddSJed Brown 
447af33a6ddSJed Brown 
448af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (fill remote requests) */
449af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
4500298fd71SBarry Smith   ierr = PetscInfo1(NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);CHKERRQ(ierr);
451af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
452af33a6ddSJed Brown     Qi = c->queueRemote[n];
453af33a6ddSJed Brown     interpIndices[0] = Qi.x;
454af33a6ddSJed Brown     interpIndices[1] = Qi.y;
455af33a6ddSJed Brown     if (c->velocityInterpLocal) {
456*61ba4676SBarry Smith       ierr = c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);CHKERRQ(ierr);
457*61ba4676SBarry Smith       ierr = c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);CHKERRQ(ierr);
458af33a6ddSJed Brown     } else {
459*61ba4676SBarry Smith       ierr = c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);CHKERRQ(ierr);
460*61ba4676SBarry Smith       ierr = c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);CHKERRQ(ierr);
461af33a6ddSJed Brown     }
462af33a6ddSJed Brown     Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
463af33a6ddSJed Brown     Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
464af33a6ddSJed Brown     c->queueRemote[n] = Qi;
465af33a6ddSJed Brown   }
466af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
467af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
468af33a6ddSJed Brown   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
469af33a6ddSJed Brown   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
470af33a6ddSJed Brown   if (c->velocityInterpLocal) {
471af33a6ddSJed Brown     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
472af33a6ddSJed Brown     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
473af33a6ddSJed Brown     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
474af33a6ddSJed Brown     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
475af33a6ddSJed Brown   }
476af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
477af33a6ddSJed Brown 
478af33a6ddSJed Brown   /* -----------------------------------------------------------------------
479af33a6ddSJed Brown      PART 2, AT t-dt
480af33a6ddSJed Brown      -----------------------------------------------------------------------*/
481af33a6ddSJed Brown 
482af33a6ddSJed Brown   /* GET POSITION AT t_n (local values) */
483af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
4840298fd71SBarry Smith   ierr = PetscInfo(NULL, "Calculating position at t_{n}\n");CHKERRQ(ierr);
485af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
486af33a6ddSJed Brown     Qi   = c->queue[n];
487af33a6ddSJed Brown     Qi.x = Qi.i - Qi.x*dt;
488af33a6ddSJed Brown     Qi.y = Qi.j - Qi.y*dt;
489af33a6ddSJed Brown 
490af33a6ddSJed Brown     /* Determine whether the position at t-dt is local */
491af33a6ddSJed Brown     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
492af33a6ddSJed Brown 
493af33a6ddSJed Brown     /* Check for Periodic boundaries and move all periodic points back onto the domain */
494af33a6ddSJed Brown     ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
495af33a6ddSJed Brown 
496af33a6ddSJed Brown     c->queue[n] = Qi;
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 = CharacteristicSendCoordinatesBegin(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 (LOCAL REQUESTS) */
505af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
506af33a6ddSJed Brown   if (c->fieldInterpLocal) {
507af33a6ddSJed Brown     ierr = DMGetLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
508af33a6ddSJed Brown     ierr = DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
509af33a6ddSJed Brown     ierr = DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
510af33a6ddSJed Brown     ierr = DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
511af33a6ddSJed Brown   }
5120298fd71SBarry Smith   ierr = PetscInfo(NULL, "Calculating local field at t_{n}\n");CHKERRQ(ierr);
513af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
514af33a6ddSJed Brown     if (c->neighbors[c->queue[n].proc] == rank) {
515af33a6ddSJed Brown       interpIndices[0] = c->queue[n].x;
516af33a6ddSJed Brown       interpIndices[1] = c->queue[n].y;
517*61ba4676SBarry Smith       if (c->fieldInterpLocal) {ierr = c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);CHKERRQ(ierr);}
518*61ba4676SBarry Smith       else {ierr = c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);CHKERRQ(ierr);}
519bbd56ea5SKarl Rupp       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
520af33a6ddSJed Brown     }
521af33a6ddSJed Brown   }
522af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
523af33a6ddSJed Brown 
524af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
525af33a6ddSJed Brown   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
526af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
527af33a6ddSJed Brown 
528af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
529af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
5300298fd71SBarry Smith   ierr = PetscInfo1(NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);CHKERRQ(ierr);
531af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
532af33a6ddSJed Brown     interpIndices[0] = c->queueRemote[n].x;
533af33a6ddSJed Brown     interpIndices[1] = c->queueRemote[n].y;
534af33a6ddSJed Brown 
535af33a6ddSJed Brown     /* for debugging purposes */
536af33a6ddSJed Brown     if (1) { /* hacked bounds test...let's do better */
537af33a6ddSJed Brown       PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1];
538af33a6ddSJed Brown 
539bbd56ea5SKarl 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);
540af33a6ddSJed Brown     }
541af33a6ddSJed Brown 
542*61ba4676SBarry Smith     if (c->fieldInterpLocal) {ierr = c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);CHKERRQ(ierr);}
543bbd56ea5SKarl Rupp     else c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
544bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
545af33a6ddSJed Brown   }
546af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
547af33a6ddSJed Brown 
548af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
549af33a6ddSJed Brown   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
550af33a6ddSJed Brown   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
551af33a6ddSJed Brown   if (c->fieldInterpLocal) {
552af33a6ddSJed Brown     ierr = DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
553af33a6ddSJed Brown     ierr = DMRestoreLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
554af33a6ddSJed Brown   }
555af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
556af33a6ddSJed Brown 
557af33a6ddSJed Brown   /* Return field of characteristics at t_n-1 */
558af33a6ddSJed Brown   ierr = PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
559af33a6ddSJed Brown   ierr = DMDAGetInfo(c->fieldDA,0,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
560af33a6ddSJed Brown   ierr = DMDAVecGetArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
561af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
562af33a6ddSJed Brown     Qi = c->queue[n];
563bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp];
564af33a6ddSJed Brown   }
565af33a6ddSJed Brown   ierr = DMDAVecRestoreArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
566af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
567af33a6ddSJed Brown   ierr = PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
568af33a6ddSJed Brown 
569af33a6ddSJed Brown   /* Cleanup */
570af33a6ddSJed Brown   ierr = PetscFree(interpIndices);CHKERRQ(ierr);
571af33a6ddSJed Brown   ierr = PetscFree(velocityValues);CHKERRQ(ierr);
572af33a6ddSJed Brown   ierr = PetscFree(velocityValuesOld);CHKERRQ(ierr);
573af33a6ddSJed Brown   ierr = PetscFree(fieldValues);CHKERRQ(ierr);
574af33a6ddSJed Brown   PetscFunctionReturn(0);
575af33a6ddSJed Brown }
576af33a6ddSJed Brown 
577af33a6ddSJed Brown #undef __FUNCT__
578af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetNeighbors"
579af33a6ddSJed Brown PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
580af33a6ddSJed Brown {
581af33a6ddSJed Brown   PetscErrorCode ierr;
582af33a6ddSJed Brown 
583af33a6ddSJed Brown   PetscFunctionBegin;
584af33a6ddSJed Brown   c->numNeighbors = numNeighbors;
5856eaa7dc8SBarry Smith   ierr = PetscFree(c->neighbors);CHKERRQ(ierr);
586785e854fSJed Brown   ierr = PetscMalloc1(numNeighbors, &c->neighbors);CHKERRQ(ierr);
587af33a6ddSJed Brown   ierr = PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));CHKERRQ(ierr);
588af33a6ddSJed Brown   PetscFunctionReturn(0);
589af33a6ddSJed Brown }
590af33a6ddSJed Brown 
591af33a6ddSJed Brown #undef __FUNCT__
592af33a6ddSJed Brown #define __FUNCT__ "CharacteristicAddPoint"
593af33a6ddSJed Brown PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
594af33a6ddSJed Brown {
595af33a6ddSJed Brown   PetscFunctionBegin;
596f23aa3ddSBarry Smith   if (c->queueSize >= c->queueMax) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax);
597af33a6ddSJed Brown   c->queue[c->queueSize++] = *point;
598af33a6ddSJed Brown   PetscFunctionReturn(0);
599af33a6ddSJed Brown }
600af33a6ddSJed Brown 
601af33a6ddSJed Brown #undef __FUNCT__
602af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSendCoordinatesBegin"
603af33a6ddSJed Brown int CharacteristicSendCoordinatesBegin(Characteristic c)
604af33a6ddSJed Brown {
605af33a6ddSJed Brown   PetscMPIInt    rank, tag = 121;
606af33a6ddSJed Brown   PetscInt       i, n;
607af33a6ddSJed Brown   PetscErrorCode ierr;
608af33a6ddSJed Brown 
609af33a6ddSJed Brown   PetscFunctionBegin;
610ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
6110b8f38c8SBarry Smith   ierr = CharacteristicHeapSort(c, c->queue, c->queueSize);CHKERRQ(ierr);
612af33a6ddSJed Brown   ierr = PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));CHKERRQ(ierr);
613bbd56ea5SKarl Rupp   for (i = 0;  i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
614af33a6ddSJed Brown   c->fillCount[0] = 0;
615af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
616ce94432eSBarry Smith     ierr = MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr);
617af33a6ddSJed Brown   }
618af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
619ce94432eSBarry Smith     ierr = MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
620af33a6ddSJed Brown   }
621af33a6ddSJed Brown   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
622af33a6ddSJed Brown   /* Initialize the remote queue */
623af33a6ddSJed Brown   c->queueLocalMax  = c->localOffsets[0]  = 0;
624af33a6ddSJed Brown   c->queueRemoteMax = c->remoteOffsets[0] = 0;
625af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
626af33a6ddSJed Brown     c->remoteOffsets[n] = c->queueRemoteMax;
627af33a6ddSJed Brown     c->queueRemoteMax  += c->fillCount[n];
628af33a6ddSJed Brown     c->localOffsets[n]  = c->queueLocalMax;
629af33a6ddSJed Brown     c->queueLocalMax   += c->needCount[n];
630af33a6ddSJed Brown   }
631af33a6ddSJed Brown   /* HACK BEGIN */
632bbd56ea5SKarl Rupp   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
633af33a6ddSJed Brown   c->needCount[0] = 0;
634af33a6ddSJed Brown   /* HACK END */
635af33a6ddSJed Brown   if (c->queueRemoteMax) {
636af33a6ddSJed Brown     ierr = PetscMalloc(sizeof(CharacteristicPointDA2D) * c->queueRemoteMax, &c->queueRemote);CHKERRQ(ierr);
6370298fd71SBarry Smith   } else c->queueRemote = NULL;
638af33a6ddSJed Brown   c->queueRemoteSize = c->queueRemoteMax;
639af33a6ddSJed Brown 
640af33a6ddSJed Brown   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
641af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
6420298fd71SBarry Smith     ierr = PetscInfo2(NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);CHKERRQ(ierr);
643ce94432eSBarry 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);
644af33a6ddSJed Brown   }
645af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
6460298fd71SBarry Smith     ierr = PetscInfo2(NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);CHKERRQ(ierr);
647ce94432eSBarry Smith     ierr = MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
648af33a6ddSJed Brown   }
649af33a6ddSJed Brown   PetscFunctionReturn(0);
650af33a6ddSJed Brown }
651af33a6ddSJed Brown 
652af33a6ddSJed Brown #undef __FUNCT__
653af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSendCoordinatesEnd"
654af33a6ddSJed Brown PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
655af33a6ddSJed Brown {
656af33a6ddSJed Brown #if 0
657af33a6ddSJed Brown   PetscMPIInt rank;
658af33a6ddSJed Brown   PetscInt    n;
659af33a6ddSJed Brown #endif
660af33a6ddSJed Brown   PetscErrorCode ierr;
661af33a6ddSJed Brown 
662af33a6ddSJed Brown   PetscFunctionBegin;
663af33a6ddSJed Brown   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
664af33a6ddSJed Brown #if 0
665ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
666af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
667f23aa3ddSBarry 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);
668af33a6ddSJed Brown   }
669af33a6ddSJed Brown #endif
670af33a6ddSJed Brown   PetscFunctionReturn(0);
671af33a6ddSJed Brown }
672af33a6ddSJed Brown 
673af33a6ddSJed Brown #undef __FUNCT__
674af33a6ddSJed Brown #define __FUNCT__ "CharacteristicGetValuesBegin"
675af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
676af33a6ddSJed Brown {
677af33a6ddSJed Brown   PetscMPIInt    tag = 121;
678af33a6ddSJed Brown   PetscInt       n;
679af33a6ddSJed Brown   PetscErrorCode ierr;
680af33a6ddSJed Brown 
681af33a6ddSJed Brown   PetscFunctionBegin;
682af33a6ddSJed Brown   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
683af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
684ce94432eSBarry 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);
685af33a6ddSJed Brown   }
686af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
687ce94432eSBarry Smith     ierr = MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
688af33a6ddSJed Brown   }
689af33a6ddSJed Brown   PetscFunctionReturn(0);
690af33a6ddSJed Brown }
691af33a6ddSJed Brown 
692af33a6ddSJed Brown #undef __FUNCT__
693af33a6ddSJed Brown #define __FUNCT__ "CharacteristicGetValuesEnd"
694af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
695af33a6ddSJed Brown {
696af33a6ddSJed Brown   PetscErrorCode ierr;
697af33a6ddSJed Brown 
698af33a6ddSJed Brown   PetscFunctionBegin;
699af33a6ddSJed Brown   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
700af33a6ddSJed Brown   /* Free queue of requests from other procs */
701af33a6ddSJed Brown   ierr = PetscFree(c->queueRemote);CHKERRQ(ierr);
702af33a6ddSJed Brown   PetscFunctionReturn(0);
703af33a6ddSJed Brown }
704af33a6ddSJed Brown 
705af33a6ddSJed Brown /*---------------------------------------------------------------------*/
706af33a6ddSJed Brown #undef __FUNCT__
7070b8f38c8SBarry Smith #define __FUNCT__ "CharacteristicHeapSort"
708af33a6ddSJed Brown /*
709af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
710af33a6ddSJed Brown */
7110b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
712af33a6ddSJed Brown /*---------------------------------------------------------------------*/
713af33a6ddSJed Brown {
7140b8f38c8SBarry Smith   PetscErrorCode          ierr;
715af33a6ddSJed Brown   CharacteristicPointDA2D temp;
7160b8f38c8SBarry Smith   PetscInt                n;
717af33a6ddSJed Brown 
7180b8f38c8SBarry Smith   PetscFunctionBegin;
719af33a6ddSJed Brown   if (0) { /* Check the order of the queue before sorting */
7200298fd71SBarry Smith     PetscInfo(NULL, "Before Heap sort\n");
721af33a6ddSJed Brown     for (n=0; n<size; n++) {
7220298fd71SBarry Smith       ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
723af33a6ddSJed Brown     }
724af33a6ddSJed Brown   }
725af33a6ddSJed Brown 
726af33a6ddSJed Brown   /* SORTING PHASE */
7270b8f38c8SBarry Smith   for (n = (size / 2)-1; n >= 0; n--) {
7280b8f38c8SBarry Smith     ierr = CharacteristicSiftDown(c, queue, n, size-1);CHKERRQ(ierr); /* Rich had size-1 here, Matt had size*/
7290b8f38c8SBarry Smith   }
730af33a6ddSJed Brown   for (n = size-1; n >= 1; n--) {
731af33a6ddSJed Brown     temp     = queue[0];
732af33a6ddSJed Brown     queue[0] = queue[n];
733af33a6ddSJed Brown     queue[n] = temp;
7340b8f38c8SBarry Smith     ierr     = CharacteristicSiftDown(c, queue, 0, n-1);CHKERRQ(ierr);
735af33a6ddSJed Brown   }
736af33a6ddSJed Brown   if (0) { /* Check the order of the queue after sorting */
7370298fd71SBarry Smith     ierr = PetscInfo(NULL, "Avter  Heap sort\n");CHKERRQ(ierr);
738af33a6ddSJed Brown     for (n=0; n<size; n++) {
7390298fd71SBarry Smith       ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
740af33a6ddSJed Brown     }
741af33a6ddSJed Brown   }
7420b8f38c8SBarry Smith   PetscFunctionReturn(0);
743af33a6ddSJed Brown }
744af33a6ddSJed Brown 
745af33a6ddSJed Brown /*---------------------------------------------------------------------*/
746af33a6ddSJed Brown #undef __FUNCT__
7470b8f38c8SBarry Smith #define __FUNCT__ "CharacteristicSiftDown"
748af33a6ddSJed Brown /*
749af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
750af33a6ddSJed Brown */
7510b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
752af33a6ddSJed Brown /*---------------------------------------------------------------------*/
753af33a6ddSJed Brown {
754af33a6ddSJed Brown   PetscBool               done = PETSC_FALSE;
755af33a6ddSJed Brown   PetscInt                maxChild;
756af33a6ddSJed Brown   CharacteristicPointDA2D temp;
757af33a6ddSJed Brown 
758af33a6ddSJed Brown   PetscFunctionBegin;
759af33a6ddSJed Brown   while ((root*2 <= bottom) && (!done)) {
760af33a6ddSJed Brown     if (root*2 == bottom) maxChild = root * 2;
761af33a6ddSJed Brown     else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2;
762af33a6ddSJed Brown     else maxChild = root * 2 + 1;
763af33a6ddSJed Brown 
764af33a6ddSJed Brown     if (queue[root].proc < queue[maxChild].proc) {
765af33a6ddSJed Brown       temp            = queue[root];
766af33a6ddSJed Brown       queue[root]     = queue[maxChild];
767af33a6ddSJed Brown       queue[maxChild] = temp;
768af33a6ddSJed Brown       root            = maxChild;
769db4deed7SKarl Rupp     } else done = PETSC_TRUE;
770af33a6ddSJed Brown   }
771af33a6ddSJed Brown   PetscFunctionReturn(0);
772af33a6ddSJed Brown }
773af33a6ddSJed Brown 
774af33a6ddSJed Brown #undef __FUNCT__
775af33a6ddSJed Brown #define __FUNCT__ "DMDAGetNeighborsRank"
776af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
777af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
778af33a6ddSJed Brown {
779bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx, by;
780af33a6ddSJed Brown   PetscBool        IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
781af33a6ddSJed Brown   MPI_Comm         comm;
782af33a6ddSJed Brown   PetscMPIInt      rank;
783af33a6ddSJed Brown   PetscInt         **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ;
784af33a6ddSJed Brown   PetscErrorCode   ierr;
785af33a6ddSJed Brown 
786af33a6ddSJed Brown   PetscFunctionBegin;
787af33a6ddSJed Brown   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
788af33a6ddSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
78922d28d08SBarry Smith   ierr = DMDAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &bx, &by,0, 0);CHKERRQ(ierr);
790af33a6ddSJed Brown 
791bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
792bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
793af33a6ddSJed Brown 
794af33a6ddSJed Brown   neighbors[0] = rank;
795af33a6ddSJed Brown   rank = 0;
796854ce69bSBarry Smith   ierr = PetscMalloc1(PJ,&procs);CHKERRQ(ierr);
797af33a6ddSJed Brown   for (pj=0; pj<PJ; pj++) {
798854ce69bSBarry Smith     ierr = PetscMalloc1(PI,&(procs[pj]));CHKERRQ(ierr);
799af33a6ddSJed Brown     for (pi=0; pi<PI; pi++) {
800af33a6ddSJed Brown       procs[pj][pi] = rank;
801af33a6ddSJed Brown       rank++;
802af33a6ddSJed Brown     }
803af33a6ddSJed Brown   }
804af33a6ddSJed Brown 
805af33a6ddSJed Brown   pi  = neighbors[0] % PI;
806af33a6ddSJed Brown   pj  = neighbors[0] / PI;
807af33a6ddSJed Brown   pim = pi-1;  if (pim<0) pim=PI-1;
808af33a6ddSJed Brown   pip = (pi+1)%PI;
809af33a6ddSJed Brown   pjm = pj-1;  if (pjm<0) pjm=PJ-1;
810af33a6ddSJed Brown   pjp = (pj+1)%PJ;
811af33a6ddSJed Brown 
812af33a6ddSJed Brown   neighbors[1] = procs[pj] [pim];
813af33a6ddSJed Brown   neighbors[2] = procs[pjp][pim];
814af33a6ddSJed Brown   neighbors[3] = procs[pjp][pi];
815af33a6ddSJed Brown   neighbors[4] = procs[pjp][pip];
816af33a6ddSJed Brown   neighbors[5] = procs[pj] [pip];
817af33a6ddSJed Brown   neighbors[6] = procs[pjm][pip];
818af33a6ddSJed Brown   neighbors[7] = procs[pjm][pi];
819af33a6ddSJed Brown   neighbors[8] = procs[pjm][pim];
820af33a6ddSJed Brown 
821af33a6ddSJed Brown   if (!IPeriodic) {
822af33a6ddSJed Brown     if (pi==0)    neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0];
823af33a6ddSJed Brown     if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0];
824af33a6ddSJed Brown   }
825af33a6ddSJed Brown 
826af33a6ddSJed Brown   if (!JPeriodic) {
827af33a6ddSJed Brown     if (pj==0)    neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0];
828af33a6ddSJed Brown     if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0];
829af33a6ddSJed Brown   }
830af33a6ddSJed Brown 
831af33a6ddSJed Brown   for (pj = 0; pj < PJ; pj++) {
832af33a6ddSJed Brown     ierr = PetscFree(procs[pj]);CHKERRQ(ierr);
833af33a6ddSJed Brown   }
834af33a6ddSJed Brown   ierr = PetscFree(procs);CHKERRQ(ierr);
835af33a6ddSJed Brown   PetscFunctionReturn(0);
836af33a6ddSJed Brown }
837af33a6ddSJed Brown 
838af33a6ddSJed Brown #undef __FUNCT__
839af33a6ddSJed Brown #define __FUNCT__ "DMDAGetNeighborRelative"
840af33a6ddSJed Brown /*
841af33a6ddSJed Brown   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
842af33a6ddSJed Brown     2 | 3 | 4
843af33a6ddSJed Brown     __|___|__
844af33a6ddSJed Brown     1 | 0 | 5
845af33a6ddSJed Brown     __|___|__
846af33a6ddSJed Brown     8 | 7 | 6
847af33a6ddSJed Brown       |   |
848af33a6ddSJed Brown */
8491cc8b206SBarry Smith PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr)
850af33a6ddSJed Brown {
851af33a6ddSJed Brown   DMDALocalInfo  info;
8521cc8b206SBarry Smith   PetscReal      is,ie,js,je;
853af33a6ddSJed Brown   PetscErrorCode ierr;
854af33a6ddSJed Brown 
855af33a6ddSJed Brown   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
8561cc8b206SBarry Smith   is   = (PetscReal) info.xs - 0.5; ie = (PetscReal) info.xs + info.xm - 0.5;
8571cc8b206SBarry Smith   js   = (PetscReal) info.ys - 0.5; je = (PetscReal) info.ys + info.ym - 0.5;
858af33a6ddSJed Brown 
859af33a6ddSJed Brown   if (ir >= is && ir <= ie) { /* center column */
860bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 0;
861bbd56ea5SKarl Rupp     else if (jr < js)         return 7;
862bbd56ea5SKarl Rupp     else                      return 3;
863af33a6ddSJed Brown   } else if (ir < is) {     /* left column */
864bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 1;
865bbd56ea5SKarl Rupp     else if (jr < js)         return 8;
866bbd56ea5SKarl Rupp     else                      return 2;
867af33a6ddSJed Brown   } else {                  /* right column */
868bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 5;
869bbd56ea5SKarl Rupp     else if (jr < js)         return 6;
870bbd56ea5SKarl Rupp     else                      return 4;
871af33a6ddSJed Brown   }
872af33a6ddSJed Brown }
873