Actual source code: characteristic.c
2: #include "../src/semiLagrange/characteristicimpl.h" /*I "characteristic.h" I*/
4: PetscCookie CHARACTERISTIC_COOKIE;
5: PetscLogEvent CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate;
6: PetscLogEvent CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange;
7: PetscLogEvent CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange;
8: PetscTruth CharacteristicRegisterAllCalled = PETSC_FALSE;
9: /*
10: Contains the list of registered characteristic routines
11: */
12: PetscFList CharacteristicList = PETSC_NULL;
14: PetscErrorCode DAGetNeighborsRank(DA, PetscMPIInt []);
15: PetscInt DAGetNeighborRelative(DA, PassiveReal, PassiveReal);
16: PetscErrorCode DAMapToPeriodicDomain(DA, PetscScalar [] );
18: PetscErrorCode HeapSort(Characteristic, Queue, PetscInt);
19: PetscErrorCode SiftDown(Characteristic, Queue, PetscInt, PetscInt);
23: PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
24: {
25: PetscTruth iascii;
30: if (!viewer) {
31: PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);
32: }
36: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &iascii);
37: if (iascii) {
38: } else {
39: if (c->ops->view) {
40: (*c->ops->view)(c, viewer);
41: }
42: }
43: return(0);
44: }
48: PetscErrorCode CharacteristicDestroy(Characteristic c)
49: {
54: if (--((PetscObject)c)->refct > 0) return(0);
56: if (c->ops->destroy) {
57: (*c->ops->destroy)(c);
58: }
59: MPI_Type_free(&c->itemType);
60: if (c->queue) {PetscFree(c->queue);}
61: if (c->queueLocal) {PetscFree(c->queueLocal);}
62: if (c->queueRemote) {PetscFree(c->queueRemote);}
63: if (c->neighbors) {PetscFree(c->neighbors);}
64: if (c->needCount) {PetscFree(c->needCount);}
65: if (c->localOffsets) {PetscFree(c->localOffsets);}
66: if (c->fillCount) {PetscFree(c->fillCount);}
67: if (c->remoteOffsets) {PetscFree(c->remoteOffsets);}
68: if (c->request) {PetscFree(c->request);}
69: if (c->status) {PetscFree(c->status);}
70: PetscLogObjectDestroy(c);
71: PetscHeaderDestroy(c);
72: return(0);
73: }
77: PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
78: {
79: Characteristic newC;
84: *c = PETSC_NULL;
85: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
86: CharacteristicInitializePackage(PETSC_NULL);
87: #endif
89: PetscHeaderCreate(newC, _p_Characteristic, struct _CharacteristicOps, CHARACTERISTIC_COOKIE, -1, "Characteristic", comm, CharacteristicDestroy, CharacteristicView);
90: PetscLogObjectCreate(newC);
91: *c = newC;
93: newC->structured = PETSC_TRUE;
94: newC->numIds = 0;
95: newC->velocityDA = PETSC_NULL;
96: newC->velocity = PETSC_NULL;
97: newC->velocityOld = PETSC_NULL;
98: newC->numVelocityComp = 0;
99: newC->velocityComp = PETSC_NULL;
100: newC->velocityInterp = PETSC_NULL;
101: newC->velocityInterpLocal = PETSC_NULL;
102: newC->velocityCtx = PETSC_NULL;
103: newC->fieldDA = PETSC_NULL;
104: newC->field = PETSC_NULL;
105: newC->numFieldComp = 0;
106: newC->fieldComp = PETSC_NULL;
107: newC->fieldInterp = PETSC_NULL;
108: newC->fieldInterpLocal = PETSC_NULL;
109: newC->fieldCtx = PETSC_NULL;
110: newC->itemType = PETSC_NULL;
111: newC->queue = PETSC_NULL;
112: newC->queueSize = 0;
113: newC->queueMax = 0;
114: newC->queueLocal = PETSC_NULL;
115: newC->queueLocalSize = 0;
116: newC->queueLocalMax = 0;
117: newC->queueRemote = PETSC_NULL;
118: newC->queueRemoteSize = 0;
119: newC->queueRemoteMax = 0;
120: newC->numNeighbors = 0;
121: newC->neighbors = PETSC_NULL;
122: newC->needCount = PETSC_NULL;
123: newC->localOffsets = PETSC_NULL;
124: newC->fillCount = PETSC_NULL;
125: newC->remoteOffsets = PETSC_NULL;
126: newC->request = PETSC_NULL;
127: newC->status = PETSC_NULL;
128: return(0);
129: }
133: /*@C
134: CharacteristicSetType - Builds Characteristic for a particular solver.
136: Collective on Characteristic
138: Input Parameters:
139: + c - the method of characteristics context
140: - type - a known method
142: Options Database Key:
143: . -characteristic_type <method> - Sets the method; use -help for a list
144: of available methods
146: Notes:
147: See "include/characteristic.h" for available methods
149: Normally, it is best to use the CharacteristicSetFromOptions() command and
150: then set the Characteristic type from the options database rather than by using
151: this routine. Using the options database provides the user with
152: maximum flexibility in evaluating the many different Krylov methods.
153: The CharacteristicSetType() routine is provided for those situations where it
154: is necessary to set the iterative solver independently of the command
155: line or options database. This might be the case, for example, when
156: the choice of iterative solver changes during the execution of the
157: program, and the user's application is taking responsibility for
158: choosing the appropriate method. In other words, this routine is
159: not for beginners.
161: Level: intermediate
163: .keywords: Characteristic, set, method
165: .seealso: CharacteristicType
167: @*/
168: PetscErrorCode CharacteristicSetType(Characteristic c, const CharacteristicType type)
169: {
170: PetscErrorCode ierr, (*r)(Characteristic);
171: PetscTruth match;
177: PetscTypeCompare((PetscObject) c, type, &match);
178: if (match) return(0);
180: if (c->data) {
181: /* destroy the old private Characteristic context */
182: (*c->ops->destroy)(c);
183: c->data = 0;
184: }
186: PetscFListFind(CharacteristicList, ((PetscObject)c)->comm,type, (void (**)(void)) &r);
187: if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
188: c->setupcalled = 0;
189: (*r)(c);
190: PetscObjectChangeTypeName((PetscObject) c, type);
191: return(0);
192: }
196: /*@
197: CharacteristicSetUp - Sets up the internal data structures for the
198: later use of an iterative solver.
200: Collective on Characteristic
202: Input Parameter:
203: . ksp - iterative context obtained from CharacteristicCreate()
205: Level: developer
207: .keywords: Characteristic, setup
209: .seealso: CharacteristicCreate(), CharacteristicSolve(), CharacteristicDestroy()
210: @*/
211: PetscErrorCode CharacteristicSetUp(Characteristic c)
212: {
218: if (!((PetscObject)c)->type_name){
219: CharacteristicSetType(c, CHARACTERISTICDA);
220: }
222: if (c->setupcalled == 2) return(0);
224: PetscLogEventBegin(CHARACTERISTIC_SetUp,c,0,0,0);
225: if (!c->setupcalled) {
226: (*c->ops->setup)(c);
227: }
228: PetscLogEventEnd(CHARACTERISTIC_SetUp,c,0,0,0);
229: c->setupcalled = 2;
230: return(0);
231: }
235: /*@C
236: CharacteristicRegister - See CharacteristicRegisterDynamic()
238: Level: advanced
239: @*/
240: PetscErrorCode CharacteristicRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Characteristic))
241: {
243: char fullname[PETSC_MAX_PATH_LEN];
246: PetscFListConcat(path,name,fullname);
247: PetscFListAdd(&CharacteristicList,sname,fullname,(void (*)(void))function);
248: return(0);
249: }
253: PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DA da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx)
254: {
256: c->velocityDA = da;
257: c->velocity = v;
258: c->velocityOld = vOld;
259: c->numVelocityComp = numComponents;
260: c->velocityComp = components;
261: c->velocityInterp = interp;
262: c->velocityCtx = ctx;
263: return(0);
264: }
268: PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DA da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx)
269: {
271: c->velocityDA = da;
272: c->velocity = v;
273: c->velocityOld = vOld;
274: c->numVelocityComp = numComponents;
275: c->velocityComp = components;
276: c->velocityInterpLocal = interp;
277: c->velocityCtx = ctx;
278: return(0);
279: }
283: PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DA da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx)
284: {
286: #if 0
287: if (numComponents > 2) {
288: SETERRQ(PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
289: }
290: #endif
291: c->fieldDA = da;
292: c->field = v;
293: c->numFieldComp = numComponents;
294: c->fieldComp = components;
295: c->fieldInterp = interp;
296: c->fieldCtx = ctx;
297: return(0);
298: }
302: PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DA da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx)
303: {
305: #if 0
306: if (numComponents > 2) {
307: SETERRQ(PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
308: }
309: #endif
310: c->fieldDA = da;
311: c->field = v;
312: c->numFieldComp = numComponents;
313: c->fieldComp = components;
314: c->fieldInterpLocal = interp;
315: c->fieldCtx = ctx;
316: return(0);
317: }
321: PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
322: {
323: CharacteristicPointDA2D Qi;
324: DA da = c->velocityDA;
325: Vec velocityLocal, velocityLocalOld;
326: Vec fieldLocal;
327: DALocalInfo info;
328: DAPeriodicType periodic_type;
329: PetscScalar **solArray;
330: void *velocityArray;
331: void *velocityArrayOld;
332: void *fieldArray;
333: PassiveScalar *interpIndices;
334: PassiveScalar *velocityValues, *velocityValuesOld;
335: PassiveScalar *fieldValues;
336: PetscMPIInt rank;
337: PetscInt dim;
338: PetscMPIInt neighbors[9];
339: PetscInt dof;
340: PetscInt gx, gy;
341: PetscInt n, ni, nj, is, ie, js, je, qs, comp;
342: PetscErrorCode ierr;
343: PetscTruth verbose = PETSC_FALSE;
346: c->queueSize = 0;
347: MPI_Comm_rank(((PetscObject)c)->comm, &rank);
348: DAGetNeighborsRank(da, neighbors);
349: CharacteristicSetNeighbors(c, 9, neighbors);
350: CharacteristicSetUp(c);
351: /* global and local grid info */
352: DAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, &periodic_type, 0);
353: DAGetLocalInfo(da, &info);
354: ni = info.mx; nj = info.my;
355: is = info.xs; ie = info.xs+info.xm;
356: js = info.ys; je = info.ys+info.ym;
357: qs = info.xm*info.ym;
358: /* Allocation */
359: PetscMalloc(dim*sizeof(PetscScalar), &interpIndices);
360: PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValues);
361: PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValuesOld);
362: PetscMalloc(c->numFieldComp*sizeof(PetscScalar), &fieldValues);
363: PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);
365: /* -----------------------------------------------------------------------
366: PART 1, AT t-dt/2
367: -----------------------------------------------------------------------*/
368: PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);
369: /* GET POSITION AT HALF TIME IN THE PAST */
370: if (c->velocityInterpLocal) {
371: DAGetLocalVector(c->velocityDA, &velocityLocal);
372: DAGetLocalVector(c->velocityDA, &velocityLocalOld);
373: DAGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);
374: DAGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);
375: DAGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);
376: DAGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);
377: DAVecGetArray(c->velocityDA, velocityLocal, &velocityArray);
378: DAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);
379: }
380: PetscInfo(PETSC_NULL, "Calculating position at t_{n - 1/2}\n");
381: for(Qi.j = js; Qi.j < je; Qi.j++) {
382: for(Qi.i = is; Qi.i < ie; Qi.i++) {
383: interpIndices[0] = Qi.i;
384: interpIndices[1] = Qi.j;
385: if (c->velocityInterpLocal) {
386: c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
387: } else {
388: c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
389: }
390: Qi.x = Qi.i - velocityValues[0]*dt/2.0;
391: Qi.y = Qi.j - velocityValues[1]*dt/2.0;
393: /* Determine whether the position at t - dt/2 is local */
394: Qi.proc = DAGetNeighborRelative(da, Qi.x, Qi.y);
396: /* Check for Periodic boundaries and move all periodic points back onto the domain */
397: DAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));
399: if (Qi.proc && verbose) {
400: printf("[%d]Remote point (%d) at n+1/2 to neighbor %d: (i:%d, j:%d) (x:%g, y:%g)\n", rank, (int)c->queueSize+1, Qi.proc, Qi.i, Qi.j, Qi.x, Qi.y);
401: }
402: CharacteristicAddPoint(c, &Qi);
403: }
404: }
405: PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);
407: PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);
408: CharacteristicSendCoordinatesBegin(c);
409: PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);
411: PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);
412: /* Calculate velocity at t_n+1/2 (local values) */
413: PetscInfo(PETSC_NULL, "Calculating local velocities at t_{n - 1/2}\n");
414: for(n = 0; n < c->queueSize; n++) {
415: Qi = c->queue[n];
416: if (c->neighbors[Qi.proc] == rank) {
417: interpIndices[0] = Qi.x;
418: interpIndices[1] = Qi.y;
419: if (c->velocityInterpLocal) {
420: c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
421: c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
422: } else {
423: c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
424: c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
425: }
426: Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
427: Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
428: }
429: c->queue[n] = Qi;
430: }
431: PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);
433: PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);
434: CharacteristicSendCoordinatesEnd(c);
435: PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);
438: /* Calculate velocity at t_n+1/2 (fill remote requests) */
439: PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);
440: PetscInfo1(PETSC_NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);
441: for(n = 0; n < c->queueRemoteSize; n++) {
442: Qi = c->queueRemote[n];
443: interpIndices[0] = Qi.x;
444: interpIndices[1] = Qi.y;
445: if (c->velocityInterpLocal) {
446: c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
447: c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
448: } else {
449: c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
450: c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
451: }
452: Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
453: Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
454: c->queueRemote[n] = Qi;
455: }
456: PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);
457: PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);
458: CharacteristicGetValuesBegin(c);
459: CharacteristicGetValuesEnd(c);
460: if (c->velocityInterpLocal) {
461: DAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray);
462: DAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);
463: DARestoreLocalVector(c->velocityDA, &velocityLocal);
464: DARestoreLocalVector(c->velocityDA, &velocityLocalOld);
465: }
466: PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);
468: /* -----------------------------------------------------------------------
469: PART 2, AT t-dt
470: -----------------------------------------------------------------------*/
472: /* GET POSITION AT t_n (local values) */
473: PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);
474: PetscInfo(PETSC_NULL, "Calculating position at t_{n}\n");
475: for(n = 0; n < c->queueSize; n++) {
476: Qi = c->queue[n];
477: Qi.x = Qi.i - Qi.x*dt;
478: Qi.y = Qi.j - Qi.y*dt;
480: /* Determine whether the position at t-dt is local */
481: Qi.proc = DAGetNeighborRelative(da, Qi.x, Qi.y);
483: /* Check for Periodic boundaries and move all periodic points back onto the domain */
484: DAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));
486: if (Qi.proc && verbose) {
487: printf("[%d]Remote point (%d) at n to neighbor %d: (i:%d, j:%d) (x:%g, y:%g)\n", rank, (int)n, Qi.proc, Qi.i, Qi.j, Qi.x, Qi.y);
488: }
489: c->queue[n] = Qi;
490: }
491: PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);
493: PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);
494: CharacteristicSendCoordinatesBegin(c);
495: PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);
497: /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
498: PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);
499: if (c->fieldInterpLocal) {
500: DAGetLocalVector(c->fieldDA, &fieldLocal);
501: DAGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);
502: DAGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);
503: DAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);
504: }
505: PetscInfo(PETSC_NULL, "Calculating local field at t_{n}\n");
506: for(n = 0; n < c->queueSize; n++) {
507: if (c->neighbors[c->queue[n].proc] == rank) {
508: interpIndices[0] = c->queue[n].x;
509: interpIndices[1] = c->queue[n].y;
510: if (c->fieldInterpLocal) {
511: c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
512: } else {
513: c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
514: }
515: for(comp = 0; comp < c->numFieldComp; comp++) {
516: c->queue[n].field[comp] = fieldValues[comp];
517: }
518: }
519: }
520: PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);
522: PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);
523: CharacteristicSendCoordinatesEnd(c);
524: PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);
526: /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
527: PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);
528: PetscInfo1(PETSC_NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);
529: for(n = 0; n < c->queueRemoteSize; n++) {
530: if (verbose) printf("[%d]Remote point: (i:%d, j:%d) (c0:%g, c1:%g)\n", rank, c->queueRemote[n].i, c->queueRemote[n].j, c->queueRemote[n].x, c->queueRemote[n].y);
531: interpIndices[0] = c->queueRemote[n].x;
532: interpIndices[1] = c->queueRemote[n].y;
534: /* for debugging purposes */
535: if (1) { /* hacked bounds test...let's do better */
536: PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1];
538: if (( im < (PetscScalar) is - 1.) || (im > (PetscScalar) ie) || (jm < (PetscScalar) js - 1.) || (jm > (PetscScalar) je)) {
539: printf("[%d]Bounds: I (%d, %d) J (%d, %d)\n", rank, (int)is, (int)ie, (int)js, (int)je);
540: SETERRQ2(PETSC_ERR_LIB, "Nonlocal point: (%g,%g)", im, jm);
541: }
542: }
544: if (c->fieldInterpLocal) {
545: c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
546: } else {
547: c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
548: }
549: for(comp = 0; comp < c->numFieldComp; comp++) {
550: c->queueRemote[n].field[comp] = fieldValues[comp];
551: }
552: }
553: PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);
555: PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);
556: CharacteristicGetValuesBegin(c);
557: CharacteristicGetValuesEnd(c);
558: if (c->fieldInterpLocal) {
559: DAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);
560: DARestoreLocalVector(c->fieldDA, &fieldLocal);
561: }
562: PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);
564: /* Return field of characteristics at t_n-1 */
565: PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);
566: DAGetInfo(c->fieldDA, 0, 0, 0, 0, 0, 0, 0, &dof, 0, 0, 0);
567: DAVecGetArray(c->fieldDA, solution, &solArray);
568: for(n = 0; n < c->queueSize; n++) {
569: Qi = c->queue[n];
570: for(comp = 0; comp < c->numFieldComp; comp++) {
571: solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp];
572: }
573: }
574: DAVecRestoreArray(c->fieldDA, solution, &solArray);
575: PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);
576: PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);
578: /* Cleanup */
579: PetscFree(interpIndices);
580: PetscFree(velocityValues);
581: PetscFree(velocityValuesOld);
582: PetscFree(fieldValues);
583: return(0);
584: }
588: PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
589: {
593: c->numNeighbors = numNeighbors;
594: PetscMalloc(numNeighbors * sizeof(PetscMPIInt), &c->neighbors);
595: PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));
596: return(0);
597: }
601: PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
602: {
604: if (c->queueSize >= c->queueMax) {
605: SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax);
606: }
607: c->queue[c->queueSize++] = *point;
608: return(0);
609: }
613: int CharacteristicSendCoordinatesBegin(Characteristic c)
614: {
615: PetscMPIInt rank, tag = 121;
616: PetscInt i, n;
620: MPI_Comm_rank(((PetscObject)c)->comm, &rank);
621: HeapSort(c, c->queue, c->queueSize);
622: PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));
623: for(i = 0; i < c->queueSize; i++) {
624: c->needCount[c->queue[i].proc]++;
625: }
626: c->fillCount[0] = 0;
627: for(n = 1; n < c->numNeighbors; n++) {
628: MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));
629: }
630: for(n = 1; n < c->numNeighbors; n++) {
631: MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, ((PetscObject)c)->comm);
632: }
633: MPI_Waitall(c->numNeighbors-1, c->request, c->status);
634: /* Initialize the remote queue */
635: c->queueLocalMax = c->localOffsets[0] = 0;
636: c->queueRemoteMax = c->remoteOffsets[0] = 0;
637: for(n = 1; n < c->numNeighbors; n++) {
638: c->remoteOffsets[n] = c->queueRemoteMax;
639: c->queueRemoteMax += c->fillCount[n];
640: c->localOffsets[n] = c->queueLocalMax;
641: c->queueLocalMax += c->needCount[n];
642: }
643: /* HACK BEGIN */
644: for(n = 1; n < c->numNeighbors; n++) {
645: c->localOffsets[n] += c->needCount[0];
646: }
647: c->needCount[0] = 0;
648: /* HACK END */
649: if (c->queueRemoteMax) {
650: PetscMalloc(sizeof(CharacteristicPointDA2D) * c->queueRemoteMax, &c->queueRemote);
651: } else {
652: c->queueRemote = PETSC_NULL;
653: }
654: c->queueRemoteSize = c->queueRemoteMax;
656: /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
657: for(n = 1; n < c->numNeighbors; n++) {
658: PetscInfo2(PETSC_NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);
659: MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));
660: }
661: for(n = 1; n < c->numNeighbors; n++) {
662: PetscInfo2(PETSC_NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);
663: MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm);
664: }
665: return(0);
666: }
670: PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
671: {
672: #if 0
673: PetscMPIInt rank;
674: PetscInt n;
675: #endif
679: MPI_Waitall(c->numNeighbors-1, c->request, c->status);
680: #if 0
681: MPI_Comm_rank(((PetscObject)c)->comm, &rank);
682: for(n = 0; n < c->queueRemoteSize; n++) {
683: if (c->neighbors[c->queueRemote[n].proc] == rank) {
684: SETERRQ2(PETSC_ERR_PLIB, "This is fucked up, n = %d proc = %d", n, c->queueRemote[n].proc);
685: }
686: }
687: #endif
688: return(0);
689: }
693: PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
694: {
695: PetscMPIInt tag = 121;
696: PetscInt n;
700: /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
701: for(n = 1; n < c->numNeighbors; n++) {
702: MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));
703: }
704: for(n = 1; n < c->numNeighbors; n++) {
705: MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm);
706: }
707: return(0);
708: }
712: PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
713: {
717: MPI_Waitall(c->numNeighbors-1, c->request, c->status);
718: /* Free queue of requests from other procs */
719: if (c->queueRemote) {
720: PetscFree(c->queueRemote);
721: }
722: return(0);
723: }
725: /*---------------------------------------------------------------------*/
728: /*
729: Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
730: */
731: int HeapSort(Characteristic c, Queue queue, PetscInt size)
732: /*---------------------------------------------------------------------*/
733: {
734: CharacteristicPointDA2D temp;
735: int n;
736:
737: if (0) { /* Check the order of the queue before sorting */
738: PetscInfo(PETSC_NULL, "Before Heap sort\n");
739: for (n=0; n<size; n++) {
740: PetscInfo2(PETSC_NULL,"%d %d\n",n,queue[n].proc);
741: }
742: }
744: /* SORTING PHASE */
745: for (n = (size / 2)-1; n >= 0; n--)
746: SiftDown(c, queue, n, size-1); /* Rich had size-1 here, Matt had size*/
747: for (n = size-1; n >= 1; n--) {
748: temp = queue[0];
749: queue[0] = queue[n];
750: queue[n] = temp;
751: SiftDown(c, queue, 0, n-1);
752: }
753: if (0) { /* Check the order of the queue after sorting */
754: PetscInfo(PETSC_NULL, "Avter Heap sort\n");
755: for (n=0; n<size; n++) {
756: PetscInfo2(PETSC_NULL,"%d %d\n",n,queue[n].proc);
757: }
758: }
759: return 0;
760: }
762: /*---------------------------------------------------------------------*/
765: /*
766: Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
767: */
768: PetscErrorCode SiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
769: /*---------------------------------------------------------------------*/
770: {
771: PetscTruth done = PETSC_FALSE;
772: PetscInt maxChild;
773: CharacteristicPointDA2D temp;
776: while ((root*2 <= bottom) && (!done)) {
777: if (root*2 == bottom) maxChild = root * 2;
778: else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2;
779: else maxChild = root * 2 + 1;
781: if (queue[root].proc < queue[maxChild].proc) {
782: temp = queue[root];
783: queue[root] = queue[maxChild];
784: queue[maxChild] = temp;
785: root = maxChild;
786: } else
787: done = PETSC_TRUE;
788: }
789: return(0);
790: }
794: /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
795: PetscErrorCode DAGetNeighborsRank(DA da, PetscMPIInt neighbors[])
796: {
797: DAPeriodicType periodic_type;
798: PetscTruth IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
799: MPI_Comm comm;
800: PetscMPIInt rank;
801: PetscInt **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ;
805: PetscObjectGetComm((PetscObject) da, &comm);
806: MPI_Comm_rank(comm, &rank);
807: DAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &periodic_type, 0);
809: if (periodic_type==DA_XPERIODIC || periodic_type==DA_XYPERIODIC) {
810: IPeriodic = PETSC_TRUE;
811: }
812: if (periodic_type==DA_YPERIODIC || periodic_type==DA_XYPERIODIC) {
813: JPeriodic = PETSC_TRUE;
814: }
816: neighbors[0] = rank;
817: rank = 0;
818: PetscMalloc(sizeof(PetscInt*)*PJ,&procs);
819: for (pj=0;pj<PJ;pj++) {
820: PetscMalloc(sizeof(PetscInt)*PI,&(procs[pj]));
821: for (pi=0;pi<PI;pi++) {
822: procs[pj][pi] = rank;
823: rank++;
824: }
825: }
826:
827: pi = neighbors[0] % PI;
828: pj = neighbors[0] / PI;
829: pim = pi-1; if (pim<0) pim=PI-1;
830: pip = (pi+1)%PI;
831: pjm = pj-1; if (pjm<0) pjm=PJ-1;
832: pjp = (pj+1)%PJ;
834: neighbors[1] = procs[pj] [pim];
835: neighbors[2] = procs[pjp][pim];
836: neighbors[3] = procs[pjp][pi];
837: neighbors[4] = procs[pjp][pip];
838: neighbors[5] = procs[pj] [pip];
839: neighbors[6] = procs[pjm][pip];
840: neighbors[7] = procs[pjm][pi];
841: neighbors[8] = procs[pjm][pim];
843: if (!IPeriodic) {
844: if (pi==0) neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0];
845: if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0];
846: }
848: if (!JPeriodic) {
849: if (pj==0) neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0];
850: if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0];
851: }
853: for(pj = 0; pj < PJ; pj++) {
854: PetscFree(procs[pj]);
855: }
856: PetscFree(procs);
857: return(0);
858: }
862: /*
863: SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
864: 2 | 3 | 4
865: __|___|__
866: 1 | 0 | 5
867: __|___|__
868: 8 | 7 | 6
869: | |
870: */
871: PetscInt DAGetNeighborRelative(DA da, PassiveReal ir, PassiveReal jr)
872: {
873: DALocalInfo info;
874: PassiveReal is,ie,js,je;
876:
877: DAGetLocalInfo(da, &info);
878: is = (PassiveReal) info.xs - 0.5; ie = (PassiveReal) info.xs + info.xm - 0.5;
879: js = (PassiveReal) info.ys - 0.5; je = (PassiveReal) info.ys + info.ym - 0.5;
880:
881: if (ir >= is && ir <= ie) { /* center column */
882: if (jr >= js && jr <= je) {
883: return 0;
884: } else if (jr < js) {
885: return 7;
886: } else {
887: return 3;
888: }
889: } else if (ir < is) { /* left column */
890: if (jr >= js && jr <= je) {
891: return 1;
892: } else if (jr < js) {
893: return 8;
894: } else {
895: return 2;
896: }
897: } else { /* right column */
898: if (jr >= js && jr <= je) {
899: return 5;
900: } else if (jr < js) {
901: return 6;
902: } else {
903: return 4;
904: }
905: }
906: }