1af0996ceSBarry Smith #include <petsc/private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/
21e25c274SJed Brown #include <petscdmda.h>
3665c2dedSJed Brown #include <petscviewer.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
1566976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt[]);
166497c311SBarry Smith static PetscMPIInt DMDAGetNeighborRelative(DM, PetscReal, PetscReal);
17af33a6ddSJed Brown
1866976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt);
1966976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt);
20af33a6ddSJed Brown
CharacteristicView(Characteristic c,PetscViewer viewer)2166976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
22d71ae5a4SJacob Faibussowitsch {
239f196a02SMartin Diehl PetscBool isascii;
24af33a6ddSJed Brown
25af33a6ddSJed Brown PetscFunctionBegin;
26af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
2748a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
28af33a6ddSJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
29af33a6ddSJed Brown PetscCheckSameComm(c, 1, viewer, 2);
30af33a6ddSJed Brown
319f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
329f196a02SMartin Diehl if (!isascii) PetscTryTypeMethod(c, view, viewer);
333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
34af33a6ddSJed Brown }
35af33a6ddSJed Brown
3626a11704SBarry Smith /*@
3726a11704SBarry Smith CharacteristicDestroy - Destroys a `Characteristic` context created with `CharacteristicCreate()`
3826a11704SBarry Smith
3926a11704SBarry Smith Collective
4026a11704SBarry Smith
4126a11704SBarry Smith Input Parameter:
4226a11704SBarry Smith . c - the `Characteristic` context
4326a11704SBarry Smith
4426a11704SBarry Smith Level: beginner
4526a11704SBarry Smith
4626a11704SBarry Smith .seealso: `Characteristic`, `CharacteristicCreate()`
4726a11704SBarry Smith @*/
CharacteristicDestroy(Characteristic * c)48d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicDestroy(Characteristic *c)
49d71ae5a4SJacob Faibussowitsch {
50af33a6ddSJed Brown PetscFunctionBegin;
513ba16761SJacob Faibussowitsch if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
526bf464f9SBarry Smith PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
53f4f49eeaSPierre Jolivet if (--((PetscObject)*c)->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
54af33a6ddSJed Brown
55213acdd3SPierre Jolivet PetscTryTypeMethod(*c, destroy);
569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&(*c)->itemType));
579566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->queue));
589566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->queueLocal));
599566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->queueRemote));
609566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->neighbors));
619566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->needCount));
629566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->localOffsets));
639566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->fillCount));
649566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->remoteOffsets));
659566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->request));
669566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->status));
679566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(c));
683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
69af33a6ddSJed Brown }
70af33a6ddSJed Brown
7126a11704SBarry Smith /*@
7226a11704SBarry Smith CharacteristicCreate - Creates a `Characteristic` context for use with the Method of Characteristics
7326a11704SBarry Smith
7426a11704SBarry Smith Collective
7526a11704SBarry Smith
7626a11704SBarry Smith Input Parameter:
7726a11704SBarry Smith . comm - MPI communicator
7826a11704SBarry Smith
7926a11704SBarry Smith Output Parameter:
8026a11704SBarry Smith . c - the `Characteristic` context
8126a11704SBarry Smith
8226a11704SBarry Smith Level: beginner
8326a11704SBarry Smith
8426a11704SBarry Smith .seealso: `Characteristic`, `CharacteristicDestroy()`
8526a11704SBarry Smith @*/
CharacteristicCreate(MPI_Comm comm,Characteristic * c)86d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
87d71ae5a4SJacob Faibussowitsch {
88af33a6ddSJed Brown Characteristic newC;
89af33a6ddSJed Brown
90af33a6ddSJed Brown PetscFunctionBegin;
914f572ea9SToby Isaac PetscAssertPointer(c, 2);
920298fd71SBarry Smith *c = NULL;
939566063dSJacob Faibussowitsch PetscCall(CharacteristicInitializePackage());
94af33a6ddSJed Brown
959566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "Characteristic", comm, CharacteristicDestroy, CharacteristicView));
96af33a6ddSJed Brown *c = newC;
97af33a6ddSJed Brown
98af33a6ddSJed Brown newC->structured = PETSC_TRUE;
99af33a6ddSJed Brown newC->numIds = 0;
1000298fd71SBarry Smith newC->velocityDA = NULL;
1010298fd71SBarry Smith newC->velocity = NULL;
1020298fd71SBarry Smith newC->velocityOld = NULL;
103af33a6ddSJed Brown newC->numVelocityComp = 0;
1040298fd71SBarry Smith newC->velocityComp = NULL;
1050298fd71SBarry Smith newC->velocityInterp = NULL;
1060298fd71SBarry Smith newC->velocityInterpLocal = NULL;
1070298fd71SBarry Smith newC->velocityCtx = NULL;
1080298fd71SBarry Smith newC->fieldDA = NULL;
1090298fd71SBarry Smith newC->field = NULL;
110af33a6ddSJed Brown newC->numFieldComp = 0;
1110298fd71SBarry Smith newC->fieldComp = NULL;
1120298fd71SBarry Smith newC->fieldInterp = NULL;
1130298fd71SBarry Smith newC->fieldInterpLocal = NULL;
1140298fd71SBarry Smith newC->fieldCtx = NULL;
1150298fd71SBarry Smith newC->itemType = 0;
1160298fd71SBarry Smith newC->queue = NULL;
117af33a6ddSJed Brown newC->queueSize = 0;
118af33a6ddSJed Brown newC->queueMax = 0;
1190298fd71SBarry Smith newC->queueLocal = NULL;
120af33a6ddSJed Brown newC->queueLocalSize = 0;
121af33a6ddSJed Brown newC->queueLocalMax = 0;
1220298fd71SBarry Smith newC->queueRemote = NULL;
123af33a6ddSJed Brown newC->queueRemoteSize = 0;
124af33a6ddSJed Brown newC->queueRemoteMax = 0;
125af33a6ddSJed Brown newC->numNeighbors = 0;
1260298fd71SBarry Smith newC->neighbors = NULL;
1270298fd71SBarry Smith newC->needCount = NULL;
1280298fd71SBarry Smith newC->localOffsets = NULL;
1290298fd71SBarry Smith newC->fillCount = NULL;
1300298fd71SBarry Smith newC->remoteOffsets = NULL;
1310298fd71SBarry Smith newC->request = NULL;
1320298fd71SBarry Smith newC->status = NULL;
1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
134af33a6ddSJed Brown }
135af33a6ddSJed Brown
136cc4c1da9SBarry Smith /*@
137af33a6ddSJed Brown CharacteristicSetType - Builds Characteristic for a particular solver.
138af33a6ddSJed Brown
139c3339decSBarry Smith Logically Collective
140af33a6ddSJed Brown
141af33a6ddSJed Brown Input Parameters:
142af33a6ddSJed Brown + c - the method of characteristics context
143af33a6ddSJed Brown - type - a known method
144af33a6ddSJed Brown
145af33a6ddSJed Brown Options Database Key:
146af33a6ddSJed Brown . -characteristic_type <method> - Sets the method; use -help for a list
147af33a6ddSJed Brown of available methods
148af33a6ddSJed Brown
149bcf0153eSBarry Smith Level: intermediate
150bcf0153eSBarry Smith
15126a11704SBarry Smith Note:
15230a76a96SBarry Smith See "include/petsccharacteristic.h" for available methods
153af33a6ddSJed Brown
1541cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicType`
155af33a6ddSJed Brown @*/
CharacteristicSetType(Characteristic c,CharacteristicType type)156d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
157d71ae5a4SJacob Faibussowitsch {
158af33a6ddSJed Brown PetscBool match;
1595f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(Characteristic);
160af33a6ddSJed Brown
161af33a6ddSJed Brown PetscFunctionBegin;
162af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
1634f572ea9SToby Isaac PetscAssertPointer(type, 2);
164af33a6ddSJed Brown
1659566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)c, type, &match));
1663ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS);
167af33a6ddSJed Brown
168af33a6ddSJed Brown if (c->data) {
169af33a6ddSJed Brown /* destroy the old private Characteristic context */
170dbbe0bcdSBarry Smith PetscUseTypeMethod(c, destroy);
1710298fd71SBarry Smith c->ops->destroy = NULL;
1722120983fSLisandro Dalcin c->data = NULL;
173af33a6ddSJed Brown }
174af33a6ddSJed Brown
1759566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(CharacteristicList, type, &r));
1766adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)c), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
177371d2eb7SMartin Diehl c->setupcalled = PETSC_FALSE;
1789566063dSJacob Faibussowitsch PetscCall((*r)(c));
1799566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)c, type));
1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
181af33a6ddSJed Brown }
182af33a6ddSJed Brown
183af33a6ddSJed Brown /*@
184af33a6ddSJed Brown CharacteristicSetUp - Sets up the internal data structures for the
18526a11704SBarry Smith later use of a `Charactoristic` .
186af33a6ddSJed Brown
187c3339decSBarry Smith Collective
188af33a6ddSJed Brown
189af33a6ddSJed Brown Input Parameter:
19026a11704SBarry Smith . c - context obtained from CharacteristicCreate()
191af33a6ddSJed Brown
192af33a6ddSJed Brown Level: developer
193af33a6ddSJed Brown
19426a11704SBarry Smith .seealso: [](ch_ts), `Characteristic`, `CharacteristicCreate()`, `CharacteristicSolve()`, `CharacteristicDestroy()`
195af33a6ddSJed Brown @*/
CharacteristicSetUp(Characteristic c)196d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetUp(Characteristic c)
197d71ae5a4SJacob Faibussowitsch {
198af33a6ddSJed Brown PetscFunctionBegin;
199af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
200af33a6ddSJed Brown
20148a46eb9SPierre Jolivet if (!((PetscObject)c)->type_name) PetscCall(CharacteristicSetType(c, CHARACTERISTICDA));
202af33a6ddSJed Brown
203371d2eb7SMartin Diehl if (c->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
204af33a6ddSJed Brown
2059566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
206ad540459SPierre Jolivet if (!c->setupcalled) PetscUseTypeMethod(c, setup);
2079566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
208371d2eb7SMartin Diehl c->setupcalled = PETSC_TRUE;
2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
210af33a6ddSJed Brown }
211af33a6ddSJed Brown
212af33a6ddSJed Brown /*@C
21326a11704SBarry Smith CharacteristicRegister - Adds an approarch to the method of characteristics package.
2141c84c290SBarry Smith
215cc4c1da9SBarry Smith Not Collective, No Fortran Support
2161c84c290SBarry Smith
2171c84c290SBarry Smith Input Parameters:
21826a11704SBarry Smith + sname - name of a new approach
2192fe279fdSBarry Smith - function - routine to create method context
2201c84c290SBarry Smith
221bcf0153eSBarry Smith Level: advanced
222bcf0153eSBarry Smith
223b43aa488SJacob Faibussowitsch Example Usage:
2241c84c290SBarry Smith .vb
225bdf89e91SBarry Smith CharacteristicRegister("my_char", MyCharCreate);
2261c84c290SBarry Smith .ve
2271c84c290SBarry Smith
2281c84c290SBarry Smith Then, your Characteristic type can be chosen with the procedural interface via
2291c84c290SBarry Smith .vb
2301c84c290SBarry Smith CharacteristicCreate(MPI_Comm, Characteristic* &char);
2311c84c290SBarry Smith CharacteristicSetType(char,"my_char");
2321c84c290SBarry Smith .ve
2331c84c290SBarry Smith or at runtime via the option
2341c84c290SBarry Smith .vb
2351c84c290SBarry Smith -characteristic_type my_char
2361c84c290SBarry Smith .ve
2371c84c290SBarry Smith
2381c84c290SBarry Smith Notes:
23926a11704SBarry Smith `CharacteristicRegister()` may be called multiple times to add several approaches.
2401c84c290SBarry Smith
2411cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicRegisterAll()`, `CharacteristicRegisterDestroy()`
242af33a6ddSJed Brown @*/
CharacteristicRegister(const char sname[],PetscErrorCode (* function)(Characteristic))243d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicRegister(const char sname[], PetscErrorCode (*function)(Characteristic))
244d71ae5a4SJacob Faibussowitsch {
245af33a6ddSJed Brown PetscFunctionBegin;
2469566063dSJacob Faibussowitsch PetscCall(CharacteristicInitializePackage());
2479566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&CharacteristicList, sname, function));
2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
249af33a6ddSJed Brown }
250af33a6ddSJed Brown
CharacteristicSetVelocityInterpolation(Characteristic c,DM da,Vec v,Vec vOld,PetscInt numComponents,PetscInt components[],PetscErrorCode (* interp)(Vec,PetscReal[],PetscInt,PetscInt[],PetscScalar[],void *),PetscCtx ctx)251*2a8381b2SBarry Smith PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), PetscCtx ctx)
252d71ae5a4SJacob Faibussowitsch {
253af33a6ddSJed Brown PetscFunctionBegin;
254af33a6ddSJed Brown c->velocityDA = da;
255af33a6ddSJed Brown c->velocity = v;
256af33a6ddSJed Brown c->velocityOld = vOld;
257af33a6ddSJed Brown c->numVelocityComp = numComponents;
258af33a6ddSJed Brown c->velocityComp = components;
259af33a6ddSJed Brown c->velocityInterp = interp;
260af33a6ddSJed Brown c->velocityCtx = ctx;
2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
262af33a6ddSJed Brown }
263af33a6ddSJed Brown
CharacteristicSetVelocityInterpolationLocal(Characteristic c,DM da,Vec v,Vec vOld,PetscInt numComponents,PetscInt components[],PetscErrorCode (* interp)(void *,PetscReal[],PetscInt,PetscInt[],PetscScalar[],void *),PetscCtx ctx)264*2a8381b2SBarry Smith PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), PetscCtx ctx)
265d71ae5a4SJacob Faibussowitsch {
266af33a6ddSJed Brown PetscFunctionBegin;
267af33a6ddSJed Brown c->velocityDA = da;
268af33a6ddSJed Brown c->velocity = v;
269af33a6ddSJed Brown c->velocityOld = vOld;
270af33a6ddSJed Brown c->numVelocityComp = numComponents;
271af33a6ddSJed Brown c->velocityComp = components;
272af33a6ddSJed Brown c->velocityInterpLocal = interp;
273af33a6ddSJed Brown c->velocityCtx = ctx;
2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
275af33a6ddSJed Brown }
276af33a6ddSJed Brown
CharacteristicSetFieldInterpolation(Characteristic c,DM da,Vec v,PetscInt numComponents,PetscInt components[],PetscErrorCode (* interp)(Vec,PetscReal[],PetscInt,PetscInt[],PetscScalar[],void *),PetscCtx ctx)277*2a8381b2SBarry Smith PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), PetscCtx ctx)
278d71ae5a4SJacob Faibussowitsch {
279af33a6ddSJed Brown PetscFunctionBegin;
280af33a6ddSJed Brown #if 0
2813c633725SBarry Smith PetscCheck(numComponents <= 2,PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
282af33a6ddSJed Brown #endif
283af33a6ddSJed Brown c->fieldDA = da;
284af33a6ddSJed Brown c->field = v;
285af33a6ddSJed Brown c->numFieldComp = numComponents;
286af33a6ddSJed Brown c->fieldComp = components;
287af33a6ddSJed Brown c->fieldInterp = interp;
288af33a6ddSJed Brown c->fieldCtx = ctx;
2893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
290af33a6ddSJed Brown }
291af33a6ddSJed Brown
CharacteristicSetFieldInterpolationLocal(Characteristic c,DM da,Vec v,PetscInt numComponents,PetscInt components[],PetscErrorCode (* interp)(void *,PetscReal[],PetscInt,PetscInt[],PetscScalar[],void *),PetscCtx ctx)292*2a8381b2SBarry Smith PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), PetscCtx ctx)
293d71ae5a4SJacob Faibussowitsch {
294af33a6ddSJed Brown PetscFunctionBegin;
295af33a6ddSJed Brown #if 0
2963c633725SBarry Smith PetscCheck(numComponents <= 2,PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
297af33a6ddSJed Brown #endif
298af33a6ddSJed Brown c->fieldDA = da;
299af33a6ddSJed Brown c->field = v;
300af33a6ddSJed Brown c->numFieldComp = numComponents;
301af33a6ddSJed Brown c->fieldComp = components;
302af33a6ddSJed Brown c->fieldInterpLocal = interp;
303af33a6ddSJed Brown c->fieldCtx = ctx;
3043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
305af33a6ddSJed Brown }
306af33a6ddSJed Brown
30726a11704SBarry Smith /*@
30826a11704SBarry Smith CharacteristicSolve - Apply the Method of Characteristics solver
30926a11704SBarry Smith
31026a11704SBarry Smith Collective
31126a11704SBarry Smith
312d6ae5217SJose E. Roman Input Parameters:
31326a11704SBarry Smith + c - context obtained from `CharacteristicCreate()`
31426a11704SBarry Smith . dt - the time-step
31526a11704SBarry Smith - solution - vector holding the solution
31626a11704SBarry Smith
31726a11704SBarry Smith Level: developer
31826a11704SBarry Smith
31926a11704SBarry Smith .seealso: [](ch_ts), `Characteristic`, `CharacteristicCreate()`, `CharacteristicDestroy()`
32026a11704SBarry Smith @*/
CharacteristicSolve(Characteristic c,PetscReal dt,Vec solution)321d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
322d71ae5a4SJacob Faibussowitsch {
323af33a6ddSJed Brown CharacteristicPointDA2D Qi;
324af33a6ddSJed Brown DM da = c->velocityDA;
325af33a6ddSJed Brown Vec velocityLocal, velocityLocalOld;
326af33a6ddSJed Brown Vec fieldLocal;
327af33a6ddSJed Brown DMDALocalInfo info;
328af33a6ddSJed Brown PetscScalar **solArray;
329af33a6ddSJed Brown void *velocityArray;
330af33a6ddSJed Brown void *velocityArrayOld;
331af33a6ddSJed Brown void *fieldArray;
3321cc8b206SBarry Smith PetscScalar *interpIndices;
3331cc8b206SBarry Smith PetscScalar *velocityValues, *velocityValuesOld;
3341cc8b206SBarry Smith PetscScalar *fieldValues;
335af33a6ddSJed Brown PetscMPIInt rank;
336af33a6ddSJed Brown PetscInt dim;
337af33a6ddSJed Brown PetscMPIInt neighbors[9];
338af33a6ddSJed Brown PetscInt dof;
339af33a6ddSJed Brown PetscInt gx, gy;
340af33a6ddSJed Brown PetscInt n, is, ie, js, je, comp;
341af33a6ddSJed Brown
342af33a6ddSJed Brown PetscFunctionBegin;
343af33a6ddSJed Brown c->queueSize = 0;
3449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
3459566063dSJacob Faibussowitsch PetscCall(DMDAGetNeighborsRank(da, neighbors));
3469566063dSJacob Faibussowitsch PetscCall(CharacteristicSetNeighbors(c, 9, neighbors));
3479566063dSJacob Faibussowitsch PetscCall(CharacteristicSetUp(c));
348af33a6ddSJed Brown /* global and local grid info */
3499566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &gx, &gy, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3509566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
3519371c9d4SSatish Balay is = info.xs;
3529371c9d4SSatish Balay ie = info.xs + info.xm;
3539371c9d4SSatish Balay js = info.ys;
3549371c9d4SSatish Balay je = info.ys + info.ym;
355af33a6ddSJed Brown /* Allocation */
3569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim, &interpIndices));
3579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValues));
3589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValuesOld));
3599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(c->numFieldComp, &fieldValues));
3609566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
361af33a6ddSJed Brown
3621d27aa22SBarry Smith /*
363af33a6ddSJed Brown PART 1, AT t-dt/2
3641d27aa22SBarry Smith */
3659566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
366af33a6ddSJed Brown /* GET POSITION AT HALF TIME IN THE PAST */
367af33a6ddSJed Brown if (c->velocityInterpLocal) {
3689566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocal));
3699566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocalOld));
3709566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3719566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3729566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3739566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3749566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray));
3759566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
376af33a6ddSJed Brown }
3779566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n"));
378af33a6ddSJed Brown for (Qi.j = js; Qi.j < je; Qi.j++) {
379af33a6ddSJed Brown for (Qi.i = is; Qi.i < ie; Qi.i++) {
380af33a6ddSJed Brown interpIndices[0] = Qi.i;
381af33a6ddSJed Brown interpIndices[1] = Qi.j;
3829566063dSJacob Faibussowitsch if (c->velocityInterpLocal) PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3839566063dSJacob Faibussowitsch else PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
384af33a6ddSJed Brown Qi.x = Qi.i - velocityValues[0] * dt / 2.0;
385af33a6ddSJed Brown Qi.y = Qi.j - velocityValues[1] * dt / 2.0;
386af33a6ddSJed Brown
387af33a6ddSJed Brown /* Determine whether the position at t - dt/2 is local */
388af33a6ddSJed Brown Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
389af33a6ddSJed Brown
390af33a6ddSJed Brown /* Check for Periodic boundaries and move all periodic points back onto the domain */
391f4f49eeaSPierre Jolivet PetscCall(DMDAMapCoordsToPeriodicDomain(da, &Qi.x, &Qi.y));
3929566063dSJacob Faibussowitsch PetscCall(CharacteristicAddPoint(c, &Qi));
393af33a6ddSJed Brown }
394af33a6ddSJed Brown }
3959566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
396af33a6ddSJed Brown
3979566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
3989566063dSJacob Faibussowitsch PetscCall(CharacteristicSendCoordinatesBegin(c));
3999566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
400af33a6ddSJed Brown
4019566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
402af33a6ddSJed Brown /* Calculate velocity at t_n+1/2 (local values) */
4039566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n"));
404af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) {
405af33a6ddSJed Brown Qi = c->queue[n];
406af33a6ddSJed Brown if (c->neighbors[Qi.proc] == rank) {
407af33a6ddSJed Brown interpIndices[0] = Qi.x;
408af33a6ddSJed Brown interpIndices[1] = Qi.y;
409af33a6ddSJed Brown if (c->velocityInterpLocal) {
4109566063dSJacob Faibussowitsch PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4119566063dSJacob Faibussowitsch PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
412af33a6ddSJed Brown } else {
4139566063dSJacob Faibussowitsch PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4149566063dSJacob Faibussowitsch PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
415af33a6ddSJed Brown }
416af33a6ddSJed Brown Qi.x = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
417af33a6ddSJed Brown Qi.y = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
418af33a6ddSJed Brown }
419af33a6ddSJed Brown c->queue[n] = Qi;
420af33a6ddSJed Brown }
4219566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
422af33a6ddSJed Brown
4239566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
4249566063dSJacob Faibussowitsch PetscCall(CharacteristicSendCoordinatesEnd(c));
4259566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
426af33a6ddSJed Brown
427af33a6ddSJed Brown /* Calculate velocity at t_n+1/2 (fill remote requests) */
4289566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
42963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote velocities at t_{n - 1/2}\n", c->queueRemoteSize));
430af33a6ddSJed Brown for (n = 0; n < c->queueRemoteSize; n++) {
431af33a6ddSJed Brown Qi = c->queueRemote[n];
432af33a6ddSJed Brown interpIndices[0] = Qi.x;
433af33a6ddSJed Brown interpIndices[1] = Qi.y;
434af33a6ddSJed Brown if (c->velocityInterpLocal) {
4359566063dSJacob Faibussowitsch PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4369566063dSJacob Faibussowitsch PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
437af33a6ddSJed Brown } else {
4389566063dSJacob Faibussowitsch PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4399566063dSJacob Faibussowitsch PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
440af33a6ddSJed Brown }
441af33a6ddSJed Brown Qi.x = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
442af33a6ddSJed Brown Qi.y = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
443af33a6ddSJed Brown c->queueRemote[n] = Qi;
444af33a6ddSJed Brown }
4459566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
4469566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
4479566063dSJacob Faibussowitsch PetscCall(CharacteristicGetValuesBegin(c));
4489566063dSJacob Faibussowitsch PetscCall(CharacteristicGetValuesEnd(c));
449af33a6ddSJed Brown if (c->velocityInterpLocal) {
4509566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray));
4519566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
4529566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocal));
4539566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocalOld));
454af33a6ddSJed Brown }
4559566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
456af33a6ddSJed Brown
4571d27aa22SBarry Smith /*
458af33a6ddSJed Brown PART 2, AT t-dt
4591d27aa22SBarry Smith */
460af33a6ddSJed Brown
461af33a6ddSJed Brown /* GET POSITION AT t_n (local values) */
4629566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
4639566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating position at t_{n}\n"));
464af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) {
465af33a6ddSJed Brown Qi = c->queue[n];
466af33a6ddSJed Brown Qi.x = Qi.i - Qi.x * dt;
467af33a6ddSJed Brown Qi.y = Qi.j - Qi.y * dt;
468af33a6ddSJed Brown
469af33a6ddSJed Brown /* Determine whether the position at t-dt is local */
470af33a6ddSJed Brown Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
471af33a6ddSJed Brown
472af33a6ddSJed Brown /* Check for Periodic boundaries and move all periodic points back onto the domain */
473f4f49eeaSPierre Jolivet PetscCall(DMDAMapCoordsToPeriodicDomain(da, &Qi.x, &Qi.y));
474af33a6ddSJed Brown
475af33a6ddSJed Brown c->queue[n] = Qi;
476af33a6ddSJed Brown }
4779566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
478af33a6ddSJed Brown
4799566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4809566063dSJacob Faibussowitsch PetscCall(CharacteristicSendCoordinatesBegin(c));
4819566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
482af33a6ddSJed Brown
483af33a6ddSJed Brown /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
4849566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
485af33a6ddSJed Brown if (c->fieldInterpLocal) {
4869566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(c->fieldDA, &fieldLocal));
4879566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4889566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4899566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray));
490af33a6ddSJed Brown }
4919566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating local field at t_{n}\n"));
492af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) {
493af33a6ddSJed Brown if (c->neighbors[c->queue[n].proc] == rank) {
494af33a6ddSJed Brown interpIndices[0] = c->queue[n].x;
495af33a6ddSJed Brown interpIndices[1] = c->queue[n].y;
4969566063dSJacob Faibussowitsch if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4979566063dSJacob Faibussowitsch else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
498bbd56ea5SKarl Rupp for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
499af33a6ddSJed Brown }
500af33a6ddSJed Brown }
5019566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
502af33a6ddSJed Brown
5039566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
5049566063dSJacob Faibussowitsch PetscCall(CharacteristicSendCoordinatesEnd(c));
5059566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
506af33a6ddSJed Brown
507af33a6ddSJed Brown /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
5089566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
50963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote field points at t_{n}\n", c->queueRemoteSize));
510af33a6ddSJed Brown for (n = 0; n < c->queueRemoteSize; n++) {
511af33a6ddSJed Brown interpIndices[0] = c->queueRemote[n].x;
512af33a6ddSJed Brown interpIndices[1] = c->queueRemote[n].y;
513af33a6ddSJed Brown
514af33a6ddSJed Brown /* for debugging purposes */
515af33a6ddSJed Brown if (1) { /* hacked bounds test...let's do better */
5169371c9d4SSatish Balay PetscScalar im = interpIndices[0];
5179371c9d4SSatish Balay PetscScalar jm = interpIndices[1];
518af33a6ddSJed Brown
51963a3b9bcSJacob Faibussowitsch PetscCheck((im >= (PetscScalar)is - 1.) && (im <= (PetscScalar)ie) && (jm >= (PetscScalar)js - 1.) && (jm <= (PetscScalar)je), PETSC_COMM_SELF, PETSC_ERR_LIB, "Nonlocal point: (%g,%g)", (double)PetscAbsScalar(im), (double)PetscAbsScalar(jm));
520af33a6ddSJed Brown }
521af33a6ddSJed Brown
5229566063dSJacob Faibussowitsch if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
5239566063dSJacob Faibussowitsch else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
524bbd56ea5SKarl Rupp for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
525af33a6ddSJed Brown }
5269566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
527af33a6ddSJed Brown
5289566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
5299566063dSJacob Faibussowitsch PetscCall(CharacteristicGetValuesBegin(c));
5309566063dSJacob Faibussowitsch PetscCall(CharacteristicGetValuesEnd(c));
531af33a6ddSJed Brown if (c->fieldInterpLocal) {
5329566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray));
5339566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(c->fieldDA, &fieldLocal));
534af33a6ddSJed Brown }
5359566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
536af33a6ddSJed Brown
537af33a6ddSJed Brown /* Return field of characteristics at t_n-1 */
5389566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5399566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(c->fieldDA, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5409566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(c->fieldDA, solution, &solArray));
541af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) {
542af33a6ddSJed Brown Qi = c->queue[n];
543bbd56ea5SKarl Rupp for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i * dof + c->fieldComp[comp]] = Qi.field[comp];
544af33a6ddSJed Brown }
5459566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(c->fieldDA, solution, &solArray));
5469566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5479566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
548af33a6ddSJed Brown
549af33a6ddSJed Brown /* Cleanup */
5509566063dSJacob Faibussowitsch PetscCall(PetscFree(interpIndices));
5519566063dSJacob Faibussowitsch PetscCall(PetscFree(velocityValues));
5529566063dSJacob Faibussowitsch PetscCall(PetscFree(velocityValuesOld));
5539566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldValues));
5543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
555af33a6ddSJed Brown }
556af33a6ddSJed Brown
CharacteristicSetNeighbors(Characteristic c,PetscInt numNeighbors,PetscMPIInt neighbors[])557d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
558d71ae5a4SJacob Faibussowitsch {
559af33a6ddSJed Brown PetscFunctionBegin;
560835f2295SStefano Zampini PetscCall(PetscMPIIntCast(numNeighbors, &c->numNeighbors));
5619566063dSJacob Faibussowitsch PetscCall(PetscFree(c->neighbors));
5629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numNeighbors, &c->neighbors));
5639566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->neighbors, neighbors, numNeighbors));
5643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
565af33a6ddSJed Brown }
566af33a6ddSJed Brown
CharacteristicAddPoint(Characteristic c,CharacteristicPointDA2D * point)567d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
568d71ae5a4SJacob Faibussowitsch {
569af33a6ddSJed Brown PetscFunctionBegin;
57063a3b9bcSJacob Faibussowitsch PetscCheck(c->queueSize < c->queueMax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exceeded maximum queue size %" PetscInt_FMT, c->queueMax);
571af33a6ddSJed Brown c->queue[c->queueSize++] = *point;
5723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
573af33a6ddSJed Brown }
574af33a6ddSJed Brown
CharacteristicSendCoordinatesBegin(Characteristic c)5753ba16761SJacob Faibussowitsch PetscErrorCode CharacteristicSendCoordinatesBegin(Characteristic c)
576d71ae5a4SJacob Faibussowitsch {
577af33a6ddSJed Brown PetscMPIInt rank, tag = 121;
578af33a6ddSJed Brown PetscInt i, n;
579af33a6ddSJed Brown
580af33a6ddSJed Brown PetscFunctionBegin;
5819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
5829566063dSJacob Faibussowitsch PetscCall(CharacteristicHeapSort(c, c->queue, c->queueSize));
5839566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(c->needCount, c->numNeighbors));
584bbd56ea5SKarl Rupp for (i = 0; i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
585af33a6ddSJed Brown c->fillCount[0] = 0;
5866497c311SBarry Smith for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPIU_Irecv(&c->fillCount[n], 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &c->request[n - 1]));
5876497c311SBarry Smith for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPIU_Send(&c->needCount[n], 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
5889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
589af33a6ddSJed Brown /* Initialize the remote queue */
590af33a6ddSJed Brown c->queueLocalMax = c->localOffsets[0] = 0;
591af33a6ddSJed Brown c->queueRemoteMax = c->remoteOffsets[0] = 0;
592af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) {
593af33a6ddSJed Brown c->remoteOffsets[n] = c->queueRemoteMax;
594af33a6ddSJed Brown c->queueRemoteMax += c->fillCount[n];
595af33a6ddSJed Brown c->localOffsets[n] = c->queueLocalMax;
596af33a6ddSJed Brown c->queueLocalMax += c->needCount[n];
597af33a6ddSJed Brown }
598af33a6ddSJed Brown /* HACK BEGIN */
599bbd56ea5SKarl Rupp for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
600af33a6ddSJed Brown c->needCount[0] = 0;
601af33a6ddSJed Brown /* HACK END */
602ac530a7eSPierre Jolivet if (c->queueRemoteMax) PetscCall(PetscMalloc1(c->queueRemoteMax, &c->queueRemote));
603ac530a7eSPierre Jolivet else c->queueRemote = NULL;
604af33a6ddSJed Brown c->queueRemoteSize = c->queueRemoteMax;
605af33a6ddSJed Brown
606af33a6ddSJed Brown /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
607af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) {
60863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Receiving %" PetscInt_FMT " requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]));
6096497c311SBarry Smith PetscCallMPI(MPIU_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &c->request[n - 1]));
610af33a6ddSJed Brown }
611af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) {
61263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Sending %" PetscInt_FMT " requests for values from proc %d\n", c->needCount[n], c->neighbors[n]));
6136497c311SBarry Smith PetscCallMPI(MPIU_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
614af33a6ddSJed Brown }
6153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
616af33a6ddSJed Brown }
617af33a6ddSJed Brown
CharacteristicSendCoordinatesEnd(Characteristic c)618d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
619d71ae5a4SJacob Faibussowitsch {
620af33a6ddSJed Brown #if 0
621af33a6ddSJed Brown PetscMPIInt rank;
622af33a6ddSJed Brown PetscInt n;
623af33a6ddSJed Brown #endif
624af33a6ddSJed Brown
625af33a6ddSJed Brown PetscFunctionBegin;
6269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
627af33a6ddSJed Brown #if 0
6289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
6293a7d0413SPierre Jolivet for (n = 0; n < c->queueRemoteSize; n++) PetscCheck(c->neighbors[c->queueRemote[n].proc] != rank,PETSC_COMM_SELF,PETSC_ERR_PLIB, "This is messed up, n = %d proc = %d", n, c->queueRemote[n].proc);
630af33a6ddSJed Brown #endif
6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
632af33a6ddSJed Brown }
633af33a6ddSJed Brown
CharacteristicGetValuesBegin(Characteristic c)634d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
635d71ae5a4SJacob Faibussowitsch {
636af33a6ddSJed Brown PetscMPIInt tag = 121;
637af33a6ddSJed Brown PetscInt n;
638af33a6ddSJed Brown
639af33a6ddSJed Brown PetscFunctionBegin;
640d5b43468SJose E. Roman /* SEND AND RECEIVE FILLED REQUESTS for velocities at t_n+1/2 */
6416497c311SBarry Smith for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPIU_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &c->request[n - 1]));
6426497c311SBarry Smith for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPIU_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
6433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
644af33a6ddSJed Brown }
645af33a6ddSJed Brown
CharacteristicGetValuesEnd(Characteristic c)646d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
647d71ae5a4SJacob Faibussowitsch {
648af33a6ddSJed Brown PetscFunctionBegin;
6499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
650af33a6ddSJed Brown /* Free queue of requests from other procs */
6519566063dSJacob Faibussowitsch PetscCall(PetscFree(c->queueRemote));
6523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
653af33a6ddSJed Brown }
654af33a6ddSJed Brown
655af33a6ddSJed Brown /*
656af33a6ddSJed Brown Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
657af33a6ddSJed Brown */
CharacteristicHeapSort(Characteristic c,Queue queue,PetscInt size)65866976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
659af33a6ddSJed Brown {
660af33a6ddSJed Brown CharacteristicPointDA2D temp;
6610b8f38c8SBarry Smith PetscInt n;
662af33a6ddSJed Brown
6630b8f38c8SBarry Smith PetscFunctionBegin;
664af33a6ddSJed Brown if (0) { /* Check the order of the queue before sorting */
6659566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Before Heap sort\n"));
66648a46eb9SPierre Jolivet for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
667af33a6ddSJed Brown }
668af33a6ddSJed Brown
669af33a6ddSJed Brown /* SORTING PHASE */
670ac530a7eSPierre Jolivet for (n = (size / 2) - 1; n >= 0; n--) PetscCall(CharacteristicSiftDown(c, queue, n, size - 1)); /* Rich had size-1 here, Matt had size*/
671af33a6ddSJed Brown for (n = size - 1; n >= 1; n--) {
672af33a6ddSJed Brown temp = queue[0];
673af33a6ddSJed Brown queue[0] = queue[n];
674af33a6ddSJed Brown queue[n] = temp;
6759566063dSJacob Faibussowitsch PetscCall(CharacteristicSiftDown(c, queue, 0, n - 1));
676af33a6ddSJed Brown }
677af33a6ddSJed Brown if (0) { /* Check the order of the queue after sorting */
6789566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Avter Heap sort\n"));
67948a46eb9SPierre Jolivet for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
680af33a6ddSJed Brown }
6813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
682af33a6ddSJed Brown }
683af33a6ddSJed Brown
684af33a6ddSJed Brown /*
685af33a6ddSJed Brown Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
686af33a6ddSJed Brown */
CharacteristicSiftDown(Characteristic c,Queue queue,PetscInt root,PetscInt bottom)68766976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
688af33a6ddSJed Brown {
689af33a6ddSJed Brown PetscBool done = PETSC_FALSE;
690af33a6ddSJed Brown PetscInt maxChild;
691af33a6ddSJed Brown CharacteristicPointDA2D temp;
692af33a6ddSJed Brown
693af33a6ddSJed Brown PetscFunctionBegin;
694af33a6ddSJed Brown while ((root * 2 <= bottom) && (!done)) {
695af33a6ddSJed Brown if (root * 2 == bottom) maxChild = root * 2;
696af33a6ddSJed Brown else if (queue[root * 2].proc > queue[root * 2 + 1].proc) maxChild = root * 2;
697af33a6ddSJed Brown else maxChild = root * 2 + 1;
698af33a6ddSJed Brown
699af33a6ddSJed Brown if (queue[root].proc < queue[maxChild].proc) {
700af33a6ddSJed Brown temp = queue[root];
701af33a6ddSJed Brown queue[root] = queue[maxChild];
702af33a6ddSJed Brown queue[maxChild] = temp;
703af33a6ddSJed Brown root = maxChild;
704db4deed7SKarl Rupp } else done = PETSC_TRUE;
705af33a6ddSJed Brown }
7063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
707af33a6ddSJed Brown }
708af33a6ddSJed Brown
709af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
DMDAGetNeighborsRank(DM da,PetscMPIInt neighbors[])71066976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
711d71ae5a4SJacob Faibussowitsch {
712bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by;
713af33a6ddSJed Brown PetscBool IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
714af33a6ddSJed Brown MPI_Comm comm;
715af33a6ddSJed Brown PetscMPIInt rank;
716835f2295SStefano Zampini PetscMPIInt **procs, pi, pj, pim, pip, pjm, pjp, PIi, PJi;
7176497c311SBarry Smith PetscInt PI, PJ;
718af33a6ddSJed Brown
719af33a6ddSJed Brown PetscFunctionBegin;
7209566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
7219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
7229566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, &PI, &PJ, NULL, NULL, NULL, &bx, &by, NULL, NULL));
723835f2295SStefano Zampini PetscCall(PetscMPIIntCast(PI, &PIi));
724835f2295SStefano Zampini PetscCall(PetscMPIIntCast(PJ, &PJi));
725bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
726bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
727af33a6ddSJed Brown
728af33a6ddSJed Brown neighbors[0] = rank;
729af33a6ddSJed Brown rank = 0;
7309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PJ, &procs));
731af33a6ddSJed Brown for (pj = 0; pj < PJ; pj++) {
732f4f49eeaSPierre Jolivet PetscCall(PetscMalloc1(PI, &procs[pj]));
733af33a6ddSJed Brown for (pi = 0; pi < PI; pi++) {
734af33a6ddSJed Brown procs[pj][pi] = rank;
735af33a6ddSJed Brown rank++;
736af33a6ddSJed Brown }
737af33a6ddSJed Brown }
738af33a6ddSJed Brown
739af33a6ddSJed Brown pi = neighbors[0] % PI;
740af33a6ddSJed Brown pj = neighbors[0] / PI;
7419371c9d4SSatish Balay pim = pi - 1;
742835f2295SStefano Zampini if (pim < 0) pim = PIi - 1;
743835f2295SStefano Zampini pip = (pi + 1) % PIi;
7449371c9d4SSatish Balay pjm = pj - 1;
745835f2295SStefano Zampini if (pjm < 0) pjm = PJi - 1;
746835f2295SStefano Zampini pjp = (pj + 1) % PJi;
747af33a6ddSJed Brown
748af33a6ddSJed Brown neighbors[1] = procs[pj][pim];
749af33a6ddSJed Brown neighbors[2] = procs[pjp][pim];
750af33a6ddSJed Brown neighbors[3] = procs[pjp][pi];
751af33a6ddSJed Brown neighbors[4] = procs[pjp][pip];
752af33a6ddSJed Brown neighbors[5] = procs[pj][pip];
753af33a6ddSJed Brown neighbors[6] = procs[pjm][pip];
754af33a6ddSJed Brown neighbors[7] = procs[pjm][pi];
755af33a6ddSJed Brown neighbors[8] = procs[pjm][pim];
756af33a6ddSJed Brown
757af33a6ddSJed Brown if (!IPeriodic) {
758af33a6ddSJed Brown if (pi == 0) neighbors[1] = neighbors[2] = neighbors[8] = neighbors[0];
759af33a6ddSJed Brown if (pi == PI - 1) neighbors[4] = neighbors[5] = neighbors[6] = neighbors[0];
760af33a6ddSJed Brown }
761af33a6ddSJed Brown
762af33a6ddSJed Brown if (!JPeriodic) {
763af33a6ddSJed Brown if (pj == 0) neighbors[6] = neighbors[7] = neighbors[8] = neighbors[0];
764af33a6ddSJed Brown if (pj == PJ - 1) neighbors[2] = neighbors[3] = neighbors[4] = neighbors[0];
765af33a6ddSJed Brown }
766af33a6ddSJed Brown
76748a46eb9SPierre Jolivet for (pj = 0; pj < PJ; pj++) PetscCall(PetscFree(procs[pj]));
7689566063dSJacob Faibussowitsch PetscCall(PetscFree(procs));
7693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
770af33a6ddSJed Brown }
771af33a6ddSJed Brown
772af33a6ddSJed Brown /*
773af33a6ddSJed Brown SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
774af33a6ddSJed Brown 2 | 3 | 4
775af33a6ddSJed Brown __|___|__
776af33a6ddSJed Brown 1 | 0 | 5
777af33a6ddSJed Brown __|___|__
778af33a6ddSJed Brown 8 | 7 | 6
779af33a6ddSJed Brown | |
780af33a6ddSJed Brown */
DMDAGetNeighborRelative(DM da,PetscReal ir,PetscReal jr)7816497c311SBarry Smith static PetscMPIInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr)
782d71ae5a4SJacob Faibussowitsch {
783af33a6ddSJed Brown DMDALocalInfo info;
7841cc8b206SBarry Smith PetscReal is, ie, js, je;
785af33a6ddSJed Brown
7863387f462SBarry Smith PetscCallAbort(PETSC_COMM_SELF, DMDAGetLocalInfo(da, &info));
7879371c9d4SSatish Balay is = (PetscReal)info.xs - 0.5;
7889371c9d4SSatish Balay ie = (PetscReal)info.xs + info.xm - 0.5;
7899371c9d4SSatish Balay js = (PetscReal)info.ys - 0.5;
7909371c9d4SSatish Balay je = (PetscReal)info.ys + info.ym - 0.5;
791af33a6ddSJed Brown
792af33a6ddSJed Brown if (ir >= is && ir <= ie) { /* center column */
793bbd56ea5SKarl Rupp if (jr >= js && jr <= je) return 0;
794bbd56ea5SKarl Rupp else if (jr < js) return 7;
795bbd56ea5SKarl Rupp else return 3;
796af33a6ddSJed Brown } else if (ir < is) { /* left column */
797bbd56ea5SKarl Rupp if (jr >= js && jr <= je) return 1;
798bbd56ea5SKarl Rupp else if (jr < js) return 8;
799bbd56ea5SKarl Rupp else return 2;
800af33a6ddSJed Brown } else { /* right column */
801bbd56ea5SKarl Rupp if (jr >= js && jr <= je) return 5;
802bbd56ea5SKarl Rupp else if (jr < js) return 6;
803bbd56ea5SKarl Rupp else return 4;
804af33a6ddSJed Brown }
805af33a6ddSJed Brown }
806