Actual source code: sagraph.c
1: #define PETSCDM_DLL
3: #include private/saimpl.h
5: typedef struct {
6: /* The following is used for mapping. */
7: SA edges;
8: PetscInt *supp; /* support -- global input indices from ownership range with nonempty images */
9: PetscInt m; /* size of suppot -- number of input indices with nonempty images */
10: PetscInt *image; /* image -- distinct global output indices */
11: PetscInt n; /* size of image -- number of distinct global output indices */
12: PetscInt *ij; /* concatenated local images of ALL local input elements (i.e., all indices from the local ownership range), sorted within each image */
13: PetscInt *ijlen; /* image segment boundaries for each local input index */
14: PetscInt maxijlen; /* max image segment size */
15: /* The following is used for binning. */
16: PetscInt *offset, *count;
17: } SAMapping_Graph;
19: /*
20: Increment ii by the number of unique elements of segment a[0,i-1] of a SORTED array a.
21: */
22: #define PetscIntArrayCountUnique(a, i, ii) \
23: { \
24: if(i) { \
25: PetscInt k = 0; \
26: ++(ii); \
27: while(++k < (i)) \
28: if ((a)[k] != (a)[k-1]) { \
29: ++ii; \
30: } \
31: } \
32: }
34: /*
35: Copy unique elements of segment a[0,i-1] of a SORTED array a, to aa[ii0,ii1-1]:
36: i is an input, and ii is an input (with value ii0) and output (ii1), counting
37: the number of unique elements copied.
38: */
39: #define PetscIntArrayCopyUnique(a, i, aa, ii)\
40: { \
41: if(i) { \
42: PetscInt k = 0; \
43: (aa)[(ii)] = (a)[k]; \
44: ++(ii); \
45: while (++k < (i)) \
46: if ((a)[k] != (a)[k-1]) { \
47: (aa)[(ii)] = (a)[k]; \
48: ++(ii); \
49: } \
50: } \
51: }
53: /*
54: Copy unique elements of segment a[0,i-1] of a SORTED array a, to aa[ii0,ii1-1]:
55: i is an input, and ii is an input (with value ii0) and output (ii1), counting
56: the number of unique elements copied. For each copied a, copy the corresponding
57: b to bb.
58: */
59: #define PetscIntArrayCopyUniqueWithScalar(a, b, i, aa, bb, ii) \
60: { \
61: if(i) { \
62: PetscInt k = 0; \
63: (aa)[(ii)] = (a)[k]; \
64: ++(ii); \
65: while (++k < (i)) \
66: if ((a)[k] != (a)[k-1]) { \
67: (aa)[(ii)] = (a)[k]; \
68: (bb)[(ii)] = (b)[k]; \
69: ++(ii); \
70: } \
71: } \
72: }
74: /*
75: Locate index i in the table table of length tablen. If i is found in table,
76: ii is its index, between 0 and tablen; otherwise, ii == -1.
77: count,last,low,high are auxiliary variables that help speed up the search when
78: it is carried out repeatedly (e.g., for all i in an array); these variables
79: should be passed back to this macro for each search iteration and not altered
80: between the macro invocations.
81: */
82: #define SAMappingGraphLocalizeIndex(tablen, table, i,ii,count,last,low,high) \
83: { \
84: PetscBool _9_found = PETSC_FALSE; \
85: /* Convert iindex to local by searching through table. */ \
86: if((count) > 0) { \
87: /* last and ii have valid previous values, that can be used to take \
88: advantage of the already known information about the table. */ \
89: if((i) > (last)) { \
90: /* lower bound is still valid, but the upper bound might not be.*/\
91: /* \
92: table is ordered, hence, is a subsequence of the integers. \
93: Thus, the distance between ind and last in table is no greater \
94: than the distance between them within the integers: ind - last. \
95: Therefore, high raised by ind-last is a valid upper bound on ind. \
96: */ \
97: (high) = PetscMin((mapg)->m, (high)+((i)-(last))); \
98: /* ii is the largest index in the table whose value does not \
99: exceed last; since i > last, i is located above ii within \
100: table */ \
101: (low) = (ii); \
102: } \
103: if((i) < (last)) { \
104: /* upper bound is still valid, but the lower bound might not be.*/ \
105: /* \
106: table is ordered, hence, is a subsequence of the integers. \
107: Thus, the distance between i and last in table is no greater \
108: than the distance between them within the integers: last - i. \
109: Therefore, low lowered by i-last is a valid upper bound on i. \
110: */ \
111: (low) = PetscMax(0,(low)+((i)-last)); \
112: /* ii is the largest index of the table entry not exceeding last; \
113: since i < last, i is located no higher than ii within table */ \
114: (high) = (ii); \
115: } \
116: }/* if((count) > 0) */ \
117: else { \
118: (low) = 0; \
119: (high) = (tablen); \
120: } \
121: (last) = (i); \
122: while((high) - (low) > 5) { \
123: (ii) = ((high)+(low))/2; \
124: if((i) < (table)[ii]) { \
125: (high) = (ii); \
126: } \
127: else { \
128: (low) = (ii); \
129: } \
130: } \
131: (ii) = (low); \
132: while((ii) < (high) && (table)[(ii)] <= (i)) { \
133: if((i) == (table)[(ii)]) { \
134: _9_found = PETSC_TRUE; \
135: break; \
136: } \
137: ++(ii); \
138: } \
139: if(!_9_found) (ii) = -1; \
140: }
141:
143: /*
144: Locate all integers from array iarr of length len in table of length tablen.
145: The indices of the located integers -- their locations in table -- are
146: stored in iiarr of length len.
147: If drop == PETSC_TRUE:
148: - if an integer is not found in table, it is omitted and upon completion
149: iilen has the number of located indices; iilen <= ilen in this case.
150: If drop == PETSC_FALE:
151: - if an integer is not found in table, it is replaced by -1; iilen == ilen
152: upon completion.
153: */
154: #define SAMappingGraphLocalizeIndices(tablen,table,ilen,iarr,iilen,iiarr,drop) \
155: do { \
156: PetscInt _10_last,_10_low = 0,_10_high = (tablen), _10_k, _10_ind; \
157: (iilen) = 0; \
158: for(_10_k = 0, _10_ind = -1, _10_last = -1; _10_k < (ilen); ++_10_k) { \
159: SAMappingGraphLocalizeIndex((tablen),(table),(iarr)[_10_k],_10_ind,_10_k,_10_last,_10_low,_10_high); \
160: if(_10_ind != -1 && !(drop)) (iiarr)[(iilen)++] = _10_ind; \
161: } \
162: } while(0)
164: #define SAMappingGraphOrderPointers_Private(i,w,j,mask,index,off,ii,ww,jj) \
165: do{ \
166: if((index) == SA_I) { \
167: (ii) = (i); \
168: if(((mask) & SA_J)&&(j)) (jj) = (j)+(off); \
169: else (jj) = PETSC_NULL; \
170: } \
171: else { \
172: (ii) = (j); \
173: if(((mask) & SA_I)&&(i)) (jj) = (i)+(off); \
174: else (jj) = PETSC_NULL; \
175: } \
176: if(((mask) & SA_W) && (w)) (ww) = (w); \
177: else (ww) = PETSC_NULL; \
178: } while(0)
184: static PetscErrorCode SAMappingGraphMap_Private(SAMapping map, PetscInt insize, const PetscInt inidxi[], const PetscScalar inval[], const PetscInt inidxj[], PetscInt *outsize, PetscInt outidxi[], PetscScalar outval[], PetscInt outidxj[], PetscInt outsizes[], PetscBool local, PetscBool drop)
185: {
186: SAMapping_Graph *mapg = (SAMapping_Graph*)map->data;
187: PetscInt i,j,k,ind=0;
190: j = 0;
191: for(i = 0; i < insize; ++i) {
192: if(!local) {
193: PetscInt last=0, low=0, high=0;
194: /* Convert to local by searching through mapg->supp. */
195: SAMappingGraphLocalizeIndex(mapg->m,mapg->supp,inidxi[i],ind,i,last,low,high);
196: }/* if(!local) */
197: else {
198: ind = inidxi[i];
199: }
200: if((ind < 0 || ind > mapg->m)){
201: if(!drop) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %D at %D not in the support", inidxi[i], i);
202: if(outsizes) outsizes[i] = 0;
203: continue;
204: }
205: if(outidxi || (inval && outval) || (inidxj && outidxj) ) {
206: for(k = mapg->ijlen[ind]; k < mapg->ijlen[ind+1]; ++k) {
207: if(outidxi) outidxi[j] = mapg->image[mapg->ij[k]];
208: if(inidxj&&outidxj) outidxj[j] = inidxj[i];
209: if(inval&&outval) outval[j] = inval[i];
210: ++j;
211: }
212: }
213: else {
214: j += mapg->ijlen[ind+1]-mapg->ijlen[ind];
215: }
216: if(outsizes) outsizes[i] = (mapg->ijlen[ind+1]-mapg->ijlen[ind]);
217: }/* for(i = 0; i < len; ++i) */
218: if(outsize) *outsize = j;
219: return(0);
220: }
225: static PetscErrorCode SAMappingGraphMapSA_Private(SAMapping map, SA inarr, SAIndex index, SA outarr, PetscInt **_outsizes, PetscBool *_own_outsizes, PetscBool local, PetscBool drop)
226: {
228: PetscInt *inidxi = PETSC_NULL, *inidxj = PETSC_NULL, outsize, *outidxi = PETSC_NULL, *outidxj = PETSC_NULL, outlen, outsizesarr[128],*outsizes;
229: PetscScalar *inval = PETSC_NULL, *outval = PETSC_NULL;
230: SALink inlink;
231: SAHunk inhunk, outhunk;
233: /**/
234: if(index != SA_I && index != SA_J)
235: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid index %D", index);
236: if(!(inarr->mask & index))
237: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "SA components %D have no index %D", inarr->mask, index);
238: /* Determine the size of the output. */
239: inlink = inarr->first;
240: outsize = 0;
241: while(inlink) {
242: inhunk = inlink->hunk;
243: SAMappingGraphOrderPointers_Private(inhunk->i,inhunk->w,inhunk->j,inhunk->mask,index,0,inidxi,inval,inidxj);
244: SAMappingGraphMap_Private(map,inhunk->length,inidxi,inval,inidxj,&outlen,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_FALSE,drop);
245: outsize += outlen;
246: inlink = inlink->next;
247: }
248: /* Get space for the output. */
249: SAHunkCreate(outsize,outarr->mask, &outhunk);
250: /* Get space for outsizes, if necessary. */
251: if(_outsizes) {
252: if(inarr->length <= 128) {
253: outsizes = outsizesarr;
254: *_own_outsizes = PETSC_FALSE;
255: }
256: else {
257: PetscMalloc(sizeof(PetscInt)*inarr->length, &outsizes);
258: *_own_outsizes = PETSC_TRUE;
259: }
260: }
261: else {
262: outsizes = PETSC_NULL;
263: }
264: /* Do the mapping. */
265: outsize = 0;
266: inlink = inarr->first;
267: while(inlink) {
268: inhunk = inlink->hunk;
269: SAMappingGraphOrderPointers_Private(inhunk->i,inhunk->w,inhunk->j,inhunk->mask,index,0,inidxi,inval,inidxj);
270: SAMappingGraphOrderPointers_Private(outhunk->i,outhunk->w,outhunk->j,outhunk->mask,index,outsize,outidxi,outval,outidxj);
271: SAMappingGraphMap_Private(map,inhunk->length,inidxi,inval,inidxj,&outlen,outidxi,outval,outidxj,outsizes,local,drop);
272: outsize += outlen;
273: if(outsizes) outsizes += inhunk->length;
274: inlink = inlink->next;
275: }
276: SAAddHunk(outarr, outhunk);
277: if(_outsizes) *_outsizes = outsizes;
278: return(0);
279: }
283: static PetscErrorCode SAMappingMap_Graph(SAMapping map, SA inarr, SAIndex index, SA outarr)
284: {
287: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
288: if(index != SA_I && index != SA_J)
289: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid index %D", index);
290: if(!(inarr->mask & index))
291: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index %D not present among SA components %D", index, inarr->mask);
292: SAMappingGraphMapSA_Private(map,inarr,index,outarr,PETSC_NULL,PETSC_NULL,PETSC_FALSE,PETSC_TRUE);
293: return(0);
294: }
298: static PetscErrorCode SAMappingMapSplit_Graph(SAMapping map, SA inarr, SAIndex index, SA *outarrs)
299: {
301: PetscInt *outsizes;
302: PetscBool own_outsizes;
303: SA outarr;
305: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
306: if(index != SA_I && index != SA_J)
307: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid index %D", index);
308: if(!(inarr->mask & index))
309: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "SA components %D don't contain index %D", inarr->mask, index);
310: SACreate(inarr->mask, &outarr);
311: SAMappingGraphMapSA_Private(map,inarr,index,outarr,&outsizes,&own_outsizes,PETSC_FALSE,PETSC_TRUE);
312: SASplit(outarr, inarr->length, outsizes, outarr->mask, outarrs);
313: SADestroy(&outarr);
314: if(own_outsizes) {
315: PetscFree(outsizes);
316: }
317: return(0);
318: }
322: static PetscErrorCode SAMappingMapLocal_Graph(SAMapping map, SA inarr, SAIndex index, SA outarr)
323: {
326: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
327: if(index != SA_I && index != SA_J)
328: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid index %D", index);
329: if(!(inarr->mask & index))
330: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index %D not present among SA components %D", index, inarr->mask);
331: SAMappingGraphMapSA_Private(map,inarr,index,outarr,PETSC_NULL,PETSC_NULL,PETSC_TRUE,PETSC_TRUE);
332: return(0);
333: }
338: static PetscErrorCode SAMappingMapSplitLocal_Graph(SAMapping map, SA inarr, SAIndex index, SA *outarrs)
339: {
341: PetscInt *outsizes;
342: PetscBool own_outsizes;
343: SA outarr;
345: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
346: if(index != SA_I && index != SA_J)
347: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid index %D", index);
348: if(!(inarr->mask & index))
349: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "SA components %D don't contain index %D", inarr->mask, index);
350: SACreate(inarr->mask, &outarr);
351: SAMappingGraphMapSA_Private(map,inarr,index,outarr,&outsizes,&own_outsizes,PETSC_TRUE,PETSC_TRUE);
352: SASplit(outarr, inarr->length, outsizes, outarr->mask, outarrs);
353: SADestroy(&outarr);
354: if(own_outsizes) {
355: PetscFree(outsizes);
356: }
357: return(0);
358: }
363: static PetscErrorCode SAMappingGraphBinSA_Private(SAMapping map, SA inarr, PetscInt index, SA outarr, const PetscInt *_outsizes[], PetscBool local, PetscBool drop)
364: {
366: SAMapping_Graph *mapg = (SAMapping_Graph*)map->data;
367: PetscInt *binoff, *bincount, i,j,k,count,ind=0, *inidxi,*inidxj=PETSC_NULL,*outidxi,*outidxj;
368: PetscScalar *inval=PETSC_NULL, *outval;
369: SALink link;
370: SAHunk inhunk, outhunk;
373: /* We'll need to count contributions to each "bin" and the offset of each bin in outidx, etc. */
374: /* Allocate the bin offset array, if necessary. */
375: if(!mapg->offset) {
376: PetscMalloc((mapg->n+1)*sizeof(PetscInt), &(mapg->offset));
377: }
378: binoff = mapg->offset;
379: if(!mapg->count) {
380: PetscMalloc(mapg->n*sizeof(PetscInt), &(mapg->count));
381: }
382: bincount = mapg->count;
383: /* Now compute bin offsets */
384: for(j = 0; j < mapg->n; ++j) {
385: binoff[j] = 0;
386: }
387: binoff[mapg->n] = 0;
388: count = 0;
389: link = inarr->first;
390: while(link) {
391: inhunk = link->hunk;
392: SAMappingGraphOrderPointers_Private(inhunk->i,inhunk->w,inhunk->j,inhunk->mask,index,0,inidxi,inval,inidxj);
393: for(i = 0; i < inhunk->length; ++i) {
394: if(!local) {
395: PetscInt low = 0, last = 0, high = 0;
396: /* Convert to local by searching through mapg->supp. */
397: SAMappingGraphLocalizeIndex(mapg->m,mapg->supp,inidxi[i],ind,count,last,low,high);
398: }/* if(!local) */
399: else {
400: ind = inidxi[i];
401: }
402: if((ind < 0 || ind > mapg->m)){
403: if(!drop) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Array index %D at %D not in the support",inidxi[i],count);
404: else goto loopend1;
405: }
406: for(k = mapg->ijlen[ind]; k < mapg->ijlen[ind+1]; ++k) {
407: ++(binoff[mapg->ij[k]+1]);
408: }
409: loopend1:
410: ++count;
411: }/* for(i = 0; i < hunk->length; ++i) */
412: for(j = 0; j < mapg->n; ++j) {
413: binoff[j+1] += binoff[j];
414: }
415: link = link->next;
416: }/* while(link) */
417: /* Allocate space for output. */
418: SAHunkCreate(binoff[mapg->n],inarr->mask,&outhunk);
419: SAMappingGraphOrderPointers_Private(outhunk->i,outhunk->w,outhunk->j,outhunk->mask,index,0,outidxi,outval,outidxj);
420: /* Now bin the input indices and values. */
421: if(outidxi || (inval && outval) || (inidxj && outidxj) ) {
422: for(j = 0; j < mapg->n; ++j) {
423: bincount[j] = 0;
424: }
425: count = 0;
426: link = inarr->first;
427: while(link) {
428: inhunk = link->hunk;
429: SAMappingGraphOrderPointers_Private(inhunk->i,inhunk->w,inhunk->j,inhunk->mask,index,0,inidxi,inval,inidxj);
430: for(i = 0; i < inhunk->length; ++i) {
431: if(!local) {
432: PetscInt low = 0, last = 0, high = 0;
433: /* Convert to local by searching through mapg->supp. */
434: SAMappingGraphLocalizeIndex(mapg->m,mapg->supp,inidxi[i],ind,count,last,low,high);
435: }/* if(!local) */
436: else {
437: ind = inidxi[i];
438: }
439: if((ind < 0 || ind > mapg->m)){
440: if(!drop) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Array index %D at %D not in the support",inidxi[i],count);
441: else goto loopend2;
442: }
443: for(k = mapg->ijlen[ind]; k < mapg->ijlen[ind+1]; ++k) {
444: j = mapg->ij[k];
445: if(outidxi) outidxi[binoff[j]+bincount[j]] = inidxi[i];
446: if(outval && inval) outval[binoff[j] +bincount[j]] = inval[i];
447: if(outidxj && inidxj) outidxj[binoff[j]+bincount[j]] = inidxj[i];
448: ++bincount[j];
449: }
450: loopend2:
451: ++count;
452: }/* for(i = 0; i < hunk->length; ++i) */
453: }/* if(outidxi || (inval && outval) || (inidxj && outidxj)) */
454: link = link->next;
455: }/* while(link) */
456: if(_outsizes) *_outsizes = bincount;
457: return(0);
458: }
462: static PetscErrorCode SAMappingBin_Graph(SAMapping map, SA inarr, SAIndex index, SA outarr)
463: {
466: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
467: if(index != SA_I && index != SA_J) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid index %D", index);
468: if(!(inarr->mask & index)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index %D not present among SA components %D", index, inarr->mask);
469: SAMappingGraphBinSA_Private(map, inarr, index, outarr, PETSC_NULL, PETSC_FALSE, PETSC_TRUE);
470: return(0);
471: }
475: static PetscErrorCode SAMappingBinSplit_Graph(SAMapping map, SA inarr, SAIndex index, SA *outarrs)
476: {
478: SAMapping_Graph *mapg = (SAMapping_Graph*)map->data;
479: const PetscInt *outsizes;
480: SA outarr;
481: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
482: if(index != SA_I && index != SA_J) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid index %D", index);
483: if(!(inarr->mask & index)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index %D not present among SA components %D", index, inarr->mask);
484: SACreate(inarr->mask, &outarr);
485: SAMappingGraphBinSA_Private(map, inarr, index, outarr, &outsizes, PETSC_FALSE, PETSC_TRUE);
486: SASplit(outarr, mapg->n, outsizes, outarr->mask, outarrs);
487: SADestroy(&outarr);
489: return(0);
490: }
494: static PetscErrorCode SAMappingBinLocal_Graph(SAMapping map, SA inarr, SAIndex index, SA outarr)
495: {
498: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
499: if(index != SA_I && index != SA_J) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid index %D", index);
500: if(!(inarr->mask & index)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index %D not present among SA components %D", index, inarr->mask);
501: SAMappingGraphBinSA_Private(map, inarr, index, outarr, PETSC_NULL, PETSC_TRUE, PETSC_TRUE);
502: return(0);
503: }
508: static PetscErrorCode SAMappingBinSplitLocal_Graph(SAMapping map, SA inarr, SAIndex index, SA *outarrs)
509: {
511: SAMapping_Graph *mapg = (SAMapping_Graph*)map->data;
512: const PetscInt *outsizes;
513: SA outarr;
515: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
516: if(index != SA_I && index != SA_J) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid index %D", index);
517: if(!(inarr->mask & index)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index %D not present among SA components %D", index, inarr->mask);
518: SACreate(inarr->mask, &outarr);
519: SAMappingGraphBinSA_Private(map, inarr, index, outarr, &outsizes, PETSC_TRUE, PETSC_TRUE);
520: SASplit(outarr, mapg->n, outsizes, outarr->mask, outarrs);
521: SADestroy(&outarr);
522: return(0);
523: }
525: /* Clears everything, but edges. */
528: static PetscErrorCode SAMappingGraphClear_Private(SAMapping map)
529: {
530: SAMapping_Graph *mapg = (SAMapping_Graph*)(map->data);
533: if(mapg->supp) {
534: PetscFree(mapg->supp);
535: }
536: if(mapg->image) {
537: PetscFree(mapg->image);
538: }
539: if(mapg->ij) {
540: PetscFree(mapg->ij);
541: }
542: if(mapg->ijlen) {
543: PetscFree(mapg->ijlen);
544: }
545: if(mapg->offset) {
546: PetscFree(mapg->offset);
547: }
548: if(mapg->count) {
549: PetscFree(mapg->count);
550: }
551: mapg->m = mapg->n = mapg->maxijlen = 0;
552: return(0);
553: }
556: /*
557: Sort ix indices, if necessary.
558: If ix duplicates exist, arrange iy indices in segments corresponding
559: to the images of the same input element. Remove iy duplicates from each
560: image segment and record the number of images in ijlen. Convert the iy
561: indices to a local numbering with the corresponding global indices stored
562: in globals.
563: */
566: static PetscErrorCode SAMappingGraphAssembleEdgesLocal_Private(SAMapping map, SA edges)
567: {
569: SAMapping_Graph *mapg = (SAMapping_Graph *) map->data;
570: PetscInt *ixidx, *iyidx, ind,start, end, i, j, totalnij, totalnij2,maxnij, nij, m,n,*ij, *ijlen, *supp, *image, len;
571: PetscBool xincreasing;
572: SA medges;
573: SAHunk hunk;
574:
576: /* Warning: in this subroutine we manipulate SA internals directly. */
577: /* Clear the old data. */
578: SAMappingGraphClear_Private(map);
579: len = edges->length;
580: if(!len) {
581: return(0);
582: }
583: /* Merge the edges, if necessary. */
584: if(edges->first != edges->last) {
585: SAMerge(edges, &medges);
586: }
587: else {
588: medges = edges;
589: }
590: if(medges->first->hunk->refcnt > 1) {
591: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpectedly high reference count on the edge SA hunk: %D", medges->first->hunk->refcnt);
592: }
593: hunk = medges->first->hunk;
594: /* We manipulate the hunk data directly. */
595: ixidx = hunk->i;
596: iyidx = hunk->j;
599: /* Determine whether ixidx is strictly increasing. */
600: ind = ixidx[0];
601: xincreasing = PETSC_TRUE;
602: for(i = 1; i < len; ++i) {
603: if(ixidx[i] <= ind) {
604: xincreasing = PETSC_FALSE;
605: break;
606: }
607: }
608: if(!xincreasing) {
609: /* ixidx is not strictly increasing, no rearrangement of ixidx is also necessary. */
610: /* sort on ixidx */
611: PetscSortIntWithArray(len,ixidx,iyidx);
612: }
613: /* Now ixidx is sorted, so we march over it and record in ijlen the number of iyidx indices corresponding to each ixidx index. */
614: m = 0;
615: totalnij = 0;
616: maxnij = 0;
617: start = 0;
618: while (start < len) {
619: end = start+1;
620: while (end < len && ixidx[end] == ixidx[start]) ++end;
621: ++(m); /* count all of ixidx[start:end-1] as a single occurence of an idx index */
622: if (end - 1 > start) { /* found 2 or more of ixidx[start] in a row */
623: /* sort the relevant portion of iy */
624: PetscSortInt(end-start,iyidx+start);
625: /* count unique elements in iyidx[start,end-1] */
626: nij = 0;
627: PetscIntArrayCountUnique(iyidx+start,end-start,nij);
628: totalnij += nij;
629: maxnij = PetscMax(maxnij, nij);
630: }
631: start = end;
632: }
633: /*
634: Now we know the size of the support -- m, and the total size of concatenated image segments -- totalnij.
635: Allocate an array for recording the support indices -- supp.
636: Allocate an array for recording the images of each support index -- ij.
637: Allocate an array for counting the number of images for each support index -- ijlen.
638: */
639: PetscMalloc(sizeof(PetscInt)*(m+1), &ijlen);
640: PetscMalloc(sizeof(PetscInt)*(totalnij),&ij);
641: PetscMalloc(sizeof(PetscInt)*m, &supp);
643: /* We now record in supp only the unique ixidx indices, and in ij the unique iyidx indices in each of the image segments. */
644: i = 0;
645: j = 0;
646: start = 0;
647: ijlen[0] = 0;
648: while (start < len) {
649: end = start+1;
650: while (end < len && ixidx[end] == ixidx[start]) ++end;
651: /* the relevant portion of iy is already sorted; copy unique iyind only. */
652: nij = 0;
653: PetscIntArrayCopyUnique(iyidx+start,end-start,ij+j,nij);
654: supp[i] = ixidx[start];
655: j += nij;
656: ++i;
657: ijlen[i] = j;
658: start = end;
659: }
660: /* (A) Construct image -- the set of unique iyidx indices. */
661: /* (B) Endow the image with a local numbering. */
662: /* (C) Then convert ij to this local numbering. */
664: /* (A) */
665: PetscSortInt(len,iyidx);
666: n = 0;
667: PetscIntArrayCountUnique(iyidx,len,n);
668: PetscMalloc(sizeof(PetscInt)*n, &image);
669: n = 0;
670: PetscIntArrayCopyUnique(iyidx,len,image,n);
671: /* (B) */
672: PetscSortInt(n, image);
673: /* (C) */
674: /*
675: This is done by going through ij and for each k in it doing a binary search in image.
676: The result of the search is ind(k) -- the position of k in image.
677: Then ind(k) replace k in ij.
678: */
679: totalnij2 = 0;
680: SAMappingGraphLocalizeIndices(n,image,totalnij,ij,totalnij2,ij,PETSC_TRUE);
681: if(totalnij!=totalnij2) {
682: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of image indices befow %D and after %D localization don't match", totalnij,totalnij2);
683: }
684: mapg->supp = supp;
685: mapg->m = m;
686: mapg->image = image;
687: mapg->n = n;
688: mapg->ij = ij;
689: mapg->ijlen = ijlen;
690: mapg->maxijlen = maxnij;
692: if(edges != medges) {
693: SADestroy(&medges);
694: }
695: return(0);
696: }
699: #undef __FUNCT__
701: PetscErrorCode SAMappingGraphAddEdgesSA(SAMapping map, SA edges)
702: {
703: SAMapping_Graph *mapg = (SAMapping_Graph*)(map->data);
705: SALink link;
706: SAHunk hunk;
709:
710: SAMappingCheckType(map, SA_MAPPING_GRAPH,1);
712: link = edges->first;
713: while(link) {
714: hunk = link->hunk;
717: link = link->next;
718: }
719: SAAddArray(mapg->edges, edges);
720: map->assembled = PETSC_FALSE;
721: return(0);
722: }
724: #undef __FUNCT__
726: PetscErrorCode SAMappingGraphAddEdgesIS(SAMapping map, IS inix, IS iniy)
727: {
729: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
730: IS ix, iy;
731: const PetscInt *ixidx, *iyidx;
732: PetscInt len;
734: SAMappingCheckType(map, SA_MAPPING_GRAPH,1);
735: /*
736: if input or output vertices are not defined, assume they are the total domain or range.
737: */
738: if(!inix) {
739: ISCreateStride(((PetscObject)map)->comm,map->xlayout->n,map->xlayout->rstart,1,&(ix));
740: }
741: else ix = inix;
742: if(!iniy) {
743: ISCreateStride(((PetscObject)map)->comm,map->ylayout->n,map->ylayout->rstart,1,&(iy));
744: }
745: else iy = iniy;
746: #if defined(PETSC_USE_DEBUG)
747: /* Consistency checks. */
748: /* Make sure the IS sizes are compatible */
749: {
750: PetscInt nix,niy;
751: ISGetLocalSize(ix,&nix);
752: ISGetLocalSize(iy,&niy);
753: if (nix != niy) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local IS sizes don't match");
756: }
757: #endif
758: ISGetLocalSize(ix,&len);
759: ISGetIndices(ix, &ixidx);
760: ISGetIndices(iy, &iyidx);
761: SAAddData(mapg->edges, len, ixidx, PETSC_NULL, iyidx);
762: ISRestoreIndices(ix,&ixidx);
763: ISRestoreIndices(iy,&iyidx);
764: if(!inix) {
765: ISDestroy(&ix);
766: }
767: if(!iniy) {
768: ISDestroy(&iy);
769: }
770: map->assembled = PETSC_FALSE;
771: return(0);
772: }
774: #undef __FUNCT__
776: PetscErrorCode SAMappingGraphAddEdges(SAMapping map, PetscInt len, const PetscInt *ii, const PetscInt *jj) {
777: SAMapping_Graph *mapg = (SAMapping_Graph*)(map->data);
781:
782: SAMappingCheckType(map, SA_MAPPING_GRAPH,1);
787: SAAddData(mapg->edges, len, ii, PETSC_NULL, jj);
788: map->assembled = PETSC_FALSE;
789: return(0);
790: }
793: #undef __FUNCT__
795: static PetscErrorCode SAMappingGraphGetEdgeCount_Private(SAMapping map, PetscInt *_len)
796: {
798: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
800: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
801: if(!map->assembled) {
802: SAGetLength(mapg->edges, _len);
803: }
804: else *_len = mapg->n;
805: return(0);
806: }
809: #undef __FUNCT__
811: static PetscErrorCode SAMappingGraphGetEdges_Private(SAMapping map, PetscInt *ix, PetscInt *iy)
812: {
814: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
815: PetscInt len,i,k;
817: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
818: if(!map->assembled) {
819: SAGetLength(mapg->edges, &len);
820: }
821: else len = mapg->n;
822: if(!map->assembled) {
823: SAGetData(mapg->edges, ix, PETSC_NULL, iy);
824: }
825: else {
826: for(i = 0; i < mapg->m; ++i) {
827: for(k = mapg->ijlen[i]; k < mapg->ijlen[i+1]; ++k) {
828: if(ix) ix[k] = mapg->supp[i];
829: if(iy) iy[k] = mapg->image[mapg->ij[k]];
830: }
831: }
832: }
833: return(0);
834: }
837: #undef __FUNCT__
839: PetscErrorCode SAMappingGraphGetEdges(SAMapping map, PetscInt *_len, PetscInt **_ix, PetscInt **_iy)
840: {
842: PetscInt len, *ix = PETSC_NULL, *iy = PETSC_NULL;
844: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
845: SAMappingGraphGetEdgeCount_Private(map, &len);
846: if(!len) {
847: if(_len) *_len = len;
848: if(_ix) *_ix = PETSC_NULL;
849: if(_iy) *_iy = PETSC_NULL;
850: return(0);
851: }
852: if(_len) *_len = len;
853: if(!_ix && !_iy) return(0);
854: if(_ix) {
855: PetscMalloc(sizeof(PetscInt)*len, _ix);
856: ix = *_ix;
857: }
858: if(_iy) {
859: PetscMalloc(sizeof(PetscInt)*len, _iy);
860: iy = *_iy;
861: }
862: SAMappingGraphGetEdges_Private(map, ix,iy);
863: return(0);
864: }
867: #undef __FUNCT__
869: PetscErrorCode SAMappingGraphGetEdgesIS(SAMapping map, IS *_ix, IS *_iy)
870: {
872: PetscInt len, *ixidx, *iyidx, **_ixidx = PETSC_NULL, **_iyidx = PETSC_NULL;
874: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
875: if(_ix) _ixidx = &ixidx;
876: if(_iy) _iyidx = &iyidx;
877: SAMappingGraphGetEdges(map, &len, _ixidx,_iyidx);
878: if(_ix) {
879: ISCreateGeneral(map->xlayout->comm, len, ixidx, PETSC_OWN_POINTER, _ix);
880: }
881: if(_iy) {
882: ISCreateGeneral(map->ylayout->comm, len, iyidx, PETSC_OWN_POINTER, _iy);
883: }
884: return(0);
885: }
887: #undef __FUNCT__
889: PetscErrorCode SAMappingGraphGetEdgesSA(SAMapping map, SA *_edges)
890: {
892: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
893: SA edges;
894: SAHunk hunk;
896: SAMappingCheckType(map, SA_MAPPING_GRAPH, 1);
897: if(!_edges) return(0);
898: SACreate(SA_I|SA_J, &edges);
899: if(!map->assembled) {
900: SAAddArray(edges, mapg->edges);
901: }
902: else {
903: SAGetHunk(edges, mapg->n, &hunk);
904: SAMappingGraphGetEdges_Private(map, hunk->i, hunk->j);
905: }
906: *_edges = edges;
907: return(0);
908: }
915: PetscErrorCode SAMappingAssemblyBegin_Graph(SAMapping map)
916: {
917: SAMapping_Graph *mapg = (SAMapping_Graph*)(map->data);
919: PetscMPIInt xsize;
922: SAMappingGraphClear_Private(map);
924: MPI_Comm_size(map->xlayout->comm, &xsize);
925: if(xsize > 1) {
926: SA aedges;
927: /* Assembled edges across the xlayout comm. */
928: SACreate(mapg->edges->mask, &aedges);
929: SAAssemble(mapg->edges, SA_I, map->xlayout, aedges);
930: /* Assemble edges locally. */
931: SAMappingGraphAssembleEdgesLocal_Private(map,aedges);
932: SADestroy(&aedges);
933: }
934: else {
935: /* Assemble edges locally. */
936: SAMappingGraphAssembleEdgesLocal_Private(map, mapg->edges);
937: }
938: SAClear(mapg->edges);
939: return(0);
940: }
945: PetscErrorCode SAMappingAssemblyEnd_Graph(SAMapping map)
946: {
948: /* Currently a noop */
949: return(0);
950: }
953: #undef __FUNCT__
955: PetscErrorCode SAMappingGetSupport_Graph(SAMapping map, PetscInt *_len, PetscInt **_supp)
956: {
958: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
960: if(_len) *_len = mapg->m;
961: if(_supp) {
962: PetscMalloc(sizeof(PetscInt)*mapg->m, _supp);
963: PetscMemcpy(*_supp, mapg->supp, sizeof(PetscInt)*mapg->m);
964: }
965: return(0);
966: }
968: #undef __FUNCT__
970: PetscErrorCode SAMappingGetSupportIS_Graph(SAMapping map, IS *supp)
971: {
973: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
975: ISCreateGeneral(map->xlayout->comm, mapg->m, mapg->supp, PETSC_COPY_VALUES, supp);
976: return(0);
977: }
979: #undef __FUNCT__
981: PetscErrorCode SAMappingGetSupportSA_Graph(SAMapping map, SA *_supp)
982: {
984: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
985: SA supp;
987: SACreate(SA_I, &supp);
988: SAAddData(supp, mapg->m, mapg->supp, PETSC_NULL, PETSC_NULL);
989: *_supp = supp;
990: return(0);
991: }
994: #undef __FUNCT__
996: PetscErrorCode SAMappingGetImage_Graph(SAMapping map, PetscInt *_len, PetscInt **_image)
997: {
999: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
1001: if(_len) *_len = mapg->n;
1002: if(_image) {
1003: PetscMalloc(sizeof(PetscInt)*mapg->n, _image);
1004: PetscMemcpy(*_image, mapg->image, sizeof(PetscInt)*mapg->n);
1005: }
1006: return(0);
1007: }
1009: #undef __FUNCT__
1011: PetscErrorCode SAMappingGetImageIS_Graph(SAMapping map, IS *image)
1012: {
1014: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
1016: ISCreateGeneral(map->ylayout->comm, mapg->n, mapg->image, PETSC_COPY_VALUES, image);
1017: return(0);
1018: }
1020: #undef __FUNCT__
1022: PetscErrorCode SAMappingGetImageSA_Graph(SAMapping map, SA *_image)
1023: {
1025: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
1026: SA image;
1028: SACreate(SA_I, &image);
1029: SAAddData(image, mapg->m, mapg->image, PETSC_NULL, PETSC_NULL);
1030: *_image = image;
1031: return(0);
1032: }
1037: PetscErrorCode SAMappingGetMaxImageSize_Graph(SAMapping map, PetscInt *maxsize)
1038: {
1039: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
1041: SAMappingCheckType(map,SA_MAPPING_GRAPH,1);
1042: *maxsize = mapg->maxijlen;
1043: return(0);
1044: }
1050: PetscErrorCode SAMappingGetOperator_Graph(SAMapping map, Mat *mat)
1051: {
1053: PetscMPIInt xsize, ysize;
1054: Vec x, y;
1055: IS ix,iy;
1056: VecScatter scatter;
1058: SAMappingCheckType(map,SA_MAPPING_GRAPH,1);
1059:
1060: MPI_Comm_size(map->xlayout->comm, &xsize);
1061: if(xsize > 1) {
1062: VecCreateMPI(((PetscObject)mat)->comm, map->xlayout->n, map->xlayout->N, &x);
1064: }
1065: else {
1066: VecCreateSeq(PETSC_COMM_SELF, map->xlayout->n, &x);
1067: }
1068: MPI_Comm_size(map->ylayout->comm, &ysize);
1069: if(ysize > 1) {
1070: VecCreateMPI(((PetscObject)mat)->comm, map->ylayout->n, map->ylayout->N, &y);
1072: }
1073: else {
1074: VecCreateSeq(PETSC_COMM_SELF, map->ylayout->n, &y);
1075: }
1076: SAMappingGraphGetEdgesIS(map, &ix, &iy);
1077: VecScatterCreate(x,ix, y,iy, &scatter);
1078: MatCreateScatter(((PetscObject)mat)->comm, scatter, mat);
1079: ISDestroy(&ix);
1080: ISDestroy(&iy);
1081: VecDestroy(&x);
1082: VecDestroy(&y);
1083: return(0);
1084: }
1087: #undef __FUNCT__
1089: PetscErrorCode SAMappingView_Graph(SAMapping map, PetscViewer v)
1090: {
1092: SAMappingCheckType(map,SA_MAPPING_GRAPH,1);
1093: /* FIX: actually implement this */
1094: return(0);
1095: }
1099: #undef __FUNCT__
1101: PetscErrorCode SAMappingInvert_Graph(SAMapping map, SAMapping *_imap)
1102: {
1103: SAMapping imap;
1105: IS ix, iy;
1107: SAMappingCreate(((PetscObject)map)->comm, &imap);
1108: SAMappingSetSizes(imap, map->xlayout->n, map->ylayout->n, map->xlayout->N, map->ylayout->N);
1109: SAMappingGraphGetEdgesIS(map, &ix, &iy);
1110: SAMappingGraphAddEdgesIS(imap, iy,ix);
1111: ISDestroy(&ix);
1112: ISDestroy(&iy);
1113: SAMappingAssemblyBegin(imap);
1114: SAMappingAssemblyEnd(imap);
1115: return(0);
1116: }
1118: #undef __FUNCT__
1120: PetscErrorCode SAMappingPushforward_Graph_Graph(SAMapping map1, SAMapping map2, SAMapping *_map3)
1121: {
1123: SAMapping map3;
1124: PetscInt nsupp1, nsupp2, nsupp3, *supp1, *supp2, *supp3, imgsize1, *imgsizes1, imgsize2, *imgsizes2, *image1, *image2, *ixidx, *iyidx;
1125: PetscInt count, i1,i2,i1low,i1high,i2low,i2high,k;
1127: SAMappingCheckType(map1,SA_MAPPING_GRAPH,1);
1128: SAMappingCheckType(map2,SA_MAPPING_GRAPH,2);
1131: /*
1132: map3 _
1133: ... |
1134: |-----| |-----| ------> |
1135: ^ ^ |
1136: | ======> | -
1137: map1 | map1 |
1138: ... | map2 _ ... |
1139: | ... | |
1140: |-----| ------> | |-----|
1141: |
1142: -
1143: */
1144: SAMappingGetSupport(map1, &nsupp1, &supp1);
1145: SAMappingGetSupport(map2, &nsupp2, &supp2);
1146: /* Avoid computing the intersection, which may be unscalable in storage. */
1147: /*
1148: Count the number of images of the intersection of supports under the "upward" (1) and "rightward" (2) maps.
1149: It is done this way: supp1 is mapped by map2 obtaining offsets2, and supp2 is mapped by map1 obtaining offsets1.
1150: */
1151: PetscMalloc2(nsupp1, PetscInt, &imgsizes2, nsupp2, PetscInt, &imgsizes1);
1152: SAMappingGraphMap_Private(map1,nsupp2,supp2,PETSC_NULL,PETSC_NULL,&imgsize1,PETSC_NULL,PETSC_NULL,PETSC_NULL,imgsizes1,PETSC_FALSE,PETSC_FALSE);
1153: SAMappingGraphMap_Private(map2,nsupp1,supp1,PETSC_NULL,PETSC_NULL,&imgsize2,PETSC_NULL,PETSC_NULL,PETSC_NULL,imgsizes2,PETSC_FALSE,PETSC_FALSE);
1154: /* Count the number of supp1 indices with nonzero images in map2 -- that's the size of the intersection. */
1155: nsupp3 = 0;
1156: for(k = 0; k < nsupp1; ++k) nsupp3 += (imgsizes1[k]>0);
1157: #if defined(PETSC_USE_DEBUG)
1158: /* Now count the number of supp2 indices with nonzero images in map1: should be the same. */
1159: {
1160: PetscInt nsupp3_2 = 0;
1161: for(k = 0; k < nsupp2; ++k) nsupp3_2 += (imgsizes2[k]>0);
1162: if(nsupp3 != nsupp3_2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Intersections supports different in map1: %D and map2: %D", nsupp3, nsupp3_2);
1163: }
1164: #endif
1165: /* Allocate indices for the intersection. */
1166: PetscMalloc(sizeof(PetscInt)*nsupp3, &supp3);
1167: nsupp3 = 0;
1168: for(k = 0; k < nsupp2; ++k) {
1169: if(imgsizes1[k]) {
1170: supp3[nsupp3] = supp2[k];
1171: ++nsupp3;
1172: }
1173: }
1174: PetscFree(supp1);
1175: PetscFree(supp2);
1176:
1177: /*
1178: Now allocate the image storage and map the supp3 to obtain the "up" (1) and "right" (2) images.
1179: Recall that imgsizes1 are allocated for supp2, and imgsizes2 for supp1.
1180: */
1181: PetscMalloc2(imgsize1,PetscInt,&image1,imgsize2,PetscInt,&image2);
1182: SAMappingGraphMap_Private(map1,nsupp3,supp3,PETSC_NULL,PETSC_NULL,&imgsize1,image1,PETSC_NULL,PETSC_NULL,imgsizes1,PETSC_FALSE,PETSC_FALSE);
1183: SAMappingGraphMap_Private(map2,nsupp3,supp3,PETSC_NULL,PETSC_NULL,&imgsize2,image2,PETSC_NULL,PETSC_NULL,imgsizes2,PETSC_FALSE,PETSC_FALSE);
1184: PetscFree(supp3);
1186: /* Count the total number of arrows to add to the pushed forward SAMapping. */
1187: count = 0;
1188: for(k = 0; k < nsupp3; ++k) {
1189: count += (imgsizes1[k])*(imgsizes2[k]);
1190: }
1191: /* Allocate storage for the composed indices. */
1192: PetscMalloc2(count, PetscInt, &ixidx, count, PetscInt, &iyidx);
1193: count= 0;
1194: i1low = 0;
1195: i2low = 0;
1196: for(k = 0; k < nsupp3; ++k) {
1197: i1high = i1low + imgsizes1[k];
1198: i2high = i2low + imgsizes2[k];
1199: for(i1 = i1low; i1 < i1high; ++i1) {
1200: for(i2 = i2low; i1 < i2high; ++i2) {
1201: ixidx[count] = image1[i1];
1202: iyidx[count] = image2[i2];
1203: ++count;
1204: }
1205: }
1206: i1low = i1high;
1207: i2low = i2high;
1208: }
1209: PetscFree2(image1,image2);
1210: PetscFree2(imgsizes1,imgsizes2);
1211: /* Now construct the new SAMapping. */
1212: SAMappingCreate(((PetscObject)map1)->comm, &map3);
1213: SAMappingSetType(map3,SA_MAPPING_GRAPH);
1214: SAMappingSetSizes(map3, map1->ylayout->n, map2->ylayout->n, map1->ylayout->N, map2->ylayout->N);
1215: SAMappingGraphAddEdges(map3,count,ixidx,iyidx);
1216: SAMappingAssemblyBegin(map3);
1217: SAMappingAssemblyEnd(map3);
1218: PetscFree2(ixidx,iyidx);
1220: *_map3 = map3;
1221: return(0);
1222: }
1228: #undef __FUNCT__
1230: PetscErrorCode SAMappingPullback_Graph_Graph(SAMapping map1, SAMapping map2, SAMapping *_map3)
1231: {
1233: SAMapping imap1,map3;
1236: SAMappingCheckType(map1,SA_MAPPING_GRAPH,1);
1237: SAMappingCheckType(map2,SA_MAPPING_GRAPH,3);
1239: /*
1241: map2 _ map2 _
1242: ... | ... |
1243: |-----| ------> | |-----| ------> |
1244: ^ | ^ |
1245: | - | -
1246: map1 | ======> map1 |
1247: ... | ... | map3 _
1248: | | ... |
1249: |-----| |-----| ------> |
1250: |
1251: -
1252: Convert this to a pushforward by inverting map1 to imap1 and then pushing map2 forward along imap1
1253: (reflect the second diagram with respect to a horizontal axis and then compare with Pushforward,
1254: of just push forward "downward".)
1256: map2 _ map2 _ map2 _
1257: ... | ... | ... |
1258: |-----| ------> | |-----| ------> | |-----| ------> |
1259: ^ | | | | |
1260: | - | - | -
1261: map1 | ======> imap1 | ======> imap1 |
1262: ... | ... | ... | map3 _
1263: | V V ... |
1264: |-----| |-----| |-----| ------> |
1265: |
1266: -
1267: */
1268: SAMappingInvert_Graph(map1, &imap1);
1269: SAMappingPushforward_Graph_Graph(imap1, map2, &map3);
1270: SAMappingDestroy(&imap1);
1271: *_map3 = map3;
1272: return(0);
1273: }
1277: #undef __FUNCT__
1279: PetscErrorCode SAMappingDestroy_Graph(SAMapping map) {
1281: SAMapping_Graph *mapg = (SAMapping_Graph *)(map->data);
1282:
1284: SAMappingGraphClear_Private(map);
1285: SADestroy(&(mapg->edges));
1286: PetscFree(mapg);
1287: map->data = PETSC_NULL;
1288:
1289: map->setup = PETSC_FALSE;
1290: map->assembled = PETSC_FALSE;
1292: PetscObjectChangeTypeName((PetscObject)map,0);
1293: PetscObjectComposeFunction((PetscObject)map,"SAMappingPullback_graph_graph_C", "",PETSC_NULL);
1294: PetscObjectComposeFunction((PetscObject)map,"SAMappingPushforward_graph_graph_C", "",PETSC_NULL);
1295: return(0);
1296: }
1299: #undef __FUNCT__
1301: PetscErrorCode SAMappingCreate_Graph(SAMapping map) {
1303: SAMapping_Graph *mapg;
1305: PetscNewLog(map, SAMapping_Graph, &mapg);
1306: SACreate(SA_I|SA_J, &(mapg->edges));
1307: map->data = (void*)mapg;
1309: map->ops->view = SAMappingView_Graph;
1310: map->ops->setup = SAMappingSetUp_SAMapping;
1311: map->ops->assemblybegin = SAMappingAssemblyBegin_Graph;
1312: map->ops->assemblyend = SAMappingAssemblyEnd_Graph;
1313: map->ops->getsupport = SAMappingGetSupport_Graph;
1314: map->ops->getsupportis = SAMappingGetSupportIS_Graph;
1315: map->ops->getsupportsa = SAMappingGetSupportSA_Graph;
1316: map->ops->getimage = SAMappingGetImage_Graph;
1317: map->ops->getimageis = SAMappingGetImageIS_Graph;
1318: map->ops->getimagesa = SAMappingGetImageSA_Graph;
1319: map->ops->getmaximagesize = SAMappingGetMaxImageSize_Graph;
1320: map->ops->maplocal = SAMappingMapLocal_Graph;
1321: map->ops->map = SAMappingMap_Graph;
1322: map->ops->binlocal = SAMappingBinLocal_Graph;
1323: map->ops->bin = SAMappingBin_Graph;
1324: map->ops->mapsplitlocal = SAMappingMapSplitLocal_Graph;
1325: map->ops->mapsplit = SAMappingMapSplit_Graph;
1326: map->ops->binsplitlocal = SAMappingBinSplitLocal_Graph;
1327: map->ops->binsplit = SAMappingBinSplit_Graph;
1328: map->ops->invert = SAMappingInvert_Graph;
1329: map->ops->getoperator = SAMappingGetOperator_Graph;
1331: map->setup = PETSC_FALSE;
1332: map->assembled = PETSC_FALSE;
1334: PetscObjectComposeFunctionDynamic((PetscObject)map,
1335: "SAMappingPullback_graph_grah_C", "SAMappingPullback_ismappinggraph_ismappinggraph",
1336: SAMappingPullback_Graph_Graph);
1337: PetscObjectComposeFunctionDynamic((PetscObject)map,
1338: "SAMappingPushforward_graph_graph_C", "SAMappingPushforward_ismappinggraph_ismappinggraph",
1339: SAMappingPushforward_Graph_Graph);
1340: PetscObjectChangeTypeName((PetscObject)map, SA_MAPPING_GRAPH);
1341: return(0);
1342: }