Actual source code: sa.c
1: #define PETSCDM_DLL
3: #include private/saimpl.h
6: PetscClassId SA_MAPPING_CLASSID;
7: PetscLogEvent SA_MAPPING_Map, SA_MAPPING_MapLocal, SA_MAPPING_MapSplit, SA_MAPPING_MapSplitLocal;
8: PetscLogEvent SA_MAPPING_Bin, SA_MAPPING_BinLocal, SA_MAPPING_BinSplit, SA_MAPPING_BinSplitLocal;
9: PetscLogEvent SA_MAPPING_AssemblyBegin, SA_MAPPING_AssemblyEnd, SA_MAPPING_Invert, SA_MAPPING_Pushforward, SA_MAPPING_Pullback;
11: PetscFList SAMappingList = PETSC_NULL;
12: PetscBool SAMappingRegisterAllCalled = PETSC_FALSE;
13: PetscBool SAMappingPackageInitialized = PETSC_FALSE;
25: /*@C
26: SAMappingFinalizePackage - This function destroys everything in the SAMapping package. It is
27: called from PetscFinalize().
29: Level: developer
31: .keywords: Petsc, destroy, package
32: .seealso: PetscFinalize()
33: @*/
34: PetscErrorCode SAMappingFinalizePackage(void)
35: {
37: SAMappingPackageInitialized = PETSC_FALSE;
38: SAMappingRegisterAllCalled = PETSC_FALSE;
39: SAMappingList = PETSC_NULL;
40: return(0);
41: }
45: /*@C
46: SAMappingInitializePackage - This function initializes everything in the SAMapping package. It is called
47: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SAMappingCreate()
48: when using static libraries.
50: Input Parameter:
51: . path - The dynamic library path, or PETSC_NULL
53: Level: developer
55: .keywords: SAMapping, initialize, package
56: .seealso: PetscInitialize()
57: @*/
58: PetscErrorCode SAMappingInitializePackage(const char path[])
59: {
60: char logList[256];
61: char *className;
62: PetscBool opt;
63: PetscErrorCode ierr;
66: if (SAMappingPackageInitialized) return(0);
68: SAMappingPackageInitialized = PETSC_TRUE;
69: /* Register Classes */
70: PetscClassIdRegister("SAMapping",&SA_MAPPING_CLASSID);
71: /* Register Constructors */
72: SAMappingRegisterAll(path);
73: /* Register Events */
74: PetscLogEventRegister("SAMappingMap", SA_MAPPING_CLASSID,&SA_MAPPING_Map);
75: PetscLogEventRegister("SAMappingMapLocal", SA_MAPPING_CLASSID,&SA_MAPPING_MapLocal);
76: PetscLogEventRegister("SAMappingBin", SA_MAPPING_CLASSID,&SA_MAPPING_Bin);
77: PetscLogEventRegister("SAMappingBinLocal", SA_MAPPING_CLASSID,&SA_MAPPING_BinLocal);
78: PetscLogEventRegister("SAMappingMap", SA_MAPPING_CLASSID,&SA_MAPPING_MapSplit);
79: PetscLogEventRegister("SAMappingMapLocal", SA_MAPPING_CLASSID,&SA_MAPPING_MapSplitLocal);
80: PetscLogEventRegister("SAMappingBin", SA_MAPPING_CLASSID,&SA_MAPPING_BinSplit);
81: PetscLogEventRegister("SAMappingBinLocal", SA_MAPPING_CLASSID,&SA_MAPPING_BinSplitLocal);
82: PetscLogEventRegister("SAMappingAssemblyBegin", SA_MAPPING_CLASSID,&SA_MAPPING_AssemblyBegin);
83: PetscLogEventRegister("SAMappingAssemblyEnd", SA_MAPPING_CLASSID,&SA_MAPPING_AssemblyEnd);
84: PetscLogEventRegister("SAMappingPushforward", SA_MAPPING_CLASSID,&SA_MAPPING_Pushforward);
85: PetscLogEventRegister("SAMappingPullback", SA_MAPPING_CLASSID,&SA_MAPPING_Pullback);
86: PetscLogEventRegister("SAMappingInvert", SA_MAPPING_CLASSID,&SA_MAPPING_Invert);
88: /* Process info exclusions */
89: PetscOptionsGetString(PETSC_NULL, "-info_exclude", logList, 256, &opt);
90: if (opt) {
91: PetscStrstr(logList, "is_mapping", &className);
92: if (className) {
93: PetscInfoDeactivateClass(SA_MAPPING_CLASSID);
94: }
95: }
96: /* Process summary exclusions */
97: PetscOptionsGetString(PETSC_NULL, "-log_summary_exclude", logList, 256, &opt);
98: if (opt) {
99: PetscStrstr(logList, "is_mapping", &className);
100: if (className) {
101: PetscLogEventDeactivateClass(SA_MAPPING_CLASSID);
102: }
103: }
104: PetscRegisterFinalize(SAMappingFinalizePackage);
105: return(0);
106: }
110: /*@C
111: SAMappingRegister - See SAMappingRegisterDynamic()
113: Level: advanced
114: @*/
115: PetscErrorCode SAMappingRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SAMapping))
116: {
118: char fullname[PETSC_MAX_PATH_LEN];
121: PetscFListConcat(path,name,fullname);
122: PetscFListAdd(&SAMappingList,sname,fullname,(void (*)(void))function);
123: return(0);
124: }
129: /*@C
130: SAMappingRegisterDestroy - Frees the list of SAMapping methods that were
131: registered by SAMappingRegisterDynamic).
133: Not Collective
135: Level: developer
137: .keywords: SAMapping, register, destroy
139: .seealso: SAMappingRegisterDynamic), SAMappingRegisterAll()
140: @*/
141: PetscErrorCode SAMappingRegisterDestroy(void)
142: {
146: PetscFListDestroy(&SAMappingList);
147: SAMappingRegisterAllCalled = PETSC_FALSE;
148: return(0);
149: }
154: /*@C
155: SAMappingRegisterAll - Registers all of the mapping constructors in the SAMapping package.
157: Not Collective
159: Level: developer
161: .keywords: SAMapping, register, all
163: .seealso: SAMappingRegisterDestroy(), SAMappingRegisterDynamic), SAMappingCreate(),
164: SAMappingSetType()
165: @*/
166: PetscErrorCode SAMappingRegisterAll(const char *path)
167: {
171: SAMappingRegisterAllCalled = PETSC_TRUE;
172: SAMappingRegisterDynamic(SA_MAPPING_GRAPH,path,"SAMappingCreate_Graph",SAMappingCreate_Graph);
173: return(0);
174: }
178: /*@C
179: SAMappingMapLocal - maps an SA with local indices from the rank's support to global indices from the rank's range.
180: Since SAMapping is in general multivalued, some local indices are mapped to multiple global indices.
181: Only selected indices (I or J) are mapped; the other indices and weights, if any, are preserved on
182: the images.
185: Not collective
187: Input Parameters:
188: + map - mapping of indices
189: . inarr - input SA
190: - index - selection of the index to map (SA_I or SA_J; PETSC_NULL is equivalent to SA_I)
193: Output Parameters:
194: . outarr - SA with the selected indices mapped
197: Level: advanced
199: Concepts: mapping^indices
201: .seealso: SAMappingGetSupport(), SAMappingGetImage(), SAMappingGetSupportSizeLocal(), SAMappingGetImageSizeLocal(),
202: SAMappingMap(), SAMappingBin(), SAMappingBinLocal(), SAMappingMapSplitLocal()
204: @*/
205: PetscErrorCode SAMappingMapLocal(SAMapping map, SA inarr, SAIndex index, SA outarr)
206: {
211: if(!index) index = SA_I;
212: SAMappingCheckAssembled(map,PETSC_TRUE,1);
213: SAMappingCheckMethod(map,map->ops->maplocal,"SAMappingMapLocal");
214: PetscLogEventBegin(SA_MAPPING_MapLocal,map,0,0,0);
215: (*map->ops->maplocal)(map,inarr,index,outarr);
216: PetscLogEventEnd(SA_MAPPING_MapLocal,map,0,0,0);
217: return(0);
218: }
222: /*@C
223: SAMappingMap - maps an SA with global indices from the rank's support to global indices from the rank's range.
224: Since SAMapping is in general multivalued, some indices are mapped to multiple global indices.
225: Only indices of the selected type (I or J) are mapped; the other indices and weights, if any, are
226: preserved on the images.
228: Not collective
230: Input Parameters:
231: + map - mapping of indices
232: . inarr - input SA of indices and weights to map
233: - index - selection of the index to map (SA_I or SA_J; PETSC_NULL is equivalent to SA_I)
236: Output Parameters:
237: . outarr - SA with the selected indices mapped
240: Level: advanced
242: Concepts: mapping^indices global
244: .seealso: SAMappingGetSupport(), SAMappingGetImage(), SAMappingGetSupportSizeLocal(), SAMappingGetImageSizeLocal(),
245: SAMappingMapLocal(), SAMappingBin(), SAMappingBinLocal(), SAMappingMapSlit()
247: @*/
248: PetscErrorCode SAMappingMap(SAMapping map, SA inarr, SAIndex index, SA outarr)
249: {
254: SAMappingCheckAssembled(map,PETSC_TRUE,1);
255: if(!index) index = SA_I;
256: SAMappingCheckMethod(map,map->ops->map,"SAMappingMap");
257: PetscLogEventBegin(SA_MAPPING_Map,map,0,0,0);
258: (*map->ops->map)(map,inarr,index,outarr);
259: PetscLogEventEnd(SA_MAPPING_Map,map,0,0,0);
260: return(0);
261: }
266: /*@C
267: SAMappingBinLocal - order local indices from the rank's support into n consecutive groups or "bins" (some possibly empty)
268: according to which of the n image indices they are mapped to on this rank. The groups are concatenated
269: and returned as a single array. See SAMappingBinSplitLocal() if separate bin output is desired.
270: Since SAMapping is potentially multivalued, the same index can appear in multiple bins.
271: The binning is done on the indices of the selected type(I or J); the other indices and weights, if any,
272: are moved to the appropriate bin together with the selected indices.
275: Not collective
277: Input Parameters:
278: + map - mapping of indices
279: . array - SA with indices to bin
280: - index - selection of the index to bin on (SA_I or SA_J; PETSC_NULL is equivalent to SA_I)
283: Output Parameters:
284: . bins - SA containing concatenated binned indices; the number of bins is the same as the result of ISGetImageSizeLocal().
286: Level: advanced
288: Concepts: binning^local indices
290: .seealso: SAMappingGetSupport(), SAMappingGetImage(), SAMappingGetSupportSizeLocal(), SAMappingGetImageSizeLocal(),
291: SAMappingBin(), SAMappingMapLocal(), SAMappingMapLocal(), SAMappingBinSplitLocal()
293: @*/
294: PetscErrorCode SAMappingBinLocal(SAMapping map, SA array, SAIndex index, SA bins)
295: {
300: SAMappingCheckAssembled(map,PETSC_TRUE, 1);
301: if(!index) index = SA_I;
302: SAMappingCheckMethod(map,map->ops->binlocal,"SAMappingBinLocal");
303: PetscLogEventBegin(SA_MAPPING_BinLocal,map,0,0,0);
304: (*map->ops->binlocal)(map,array,index,bins);
305: PetscLogEventEnd(SA_MAPPING_BinLocal,map,0,0,0);
306: return(0);
307: }
312: /*@C
313: SAMappingBin - group local indices from the rank's support into n groups or "bins" (some possibly empty)
314: according to which of the n image indices they are mapped to on this rank. The groups are
315: concatenated and returned as a single array. See SAMappingBinSplit() if separate bin output
316: is desired.
317: Since SAMapping is potentially multivalued, the same index can appear in multiple bins.
318: The binning is done only on the indices of the selected type (I or J); the other indices and weights,
319: if any, are moved to the appropriate bin together with the selected indices.
322: Not collective
324: Input Parameters:
325: + map - mapping of indices
326: . array - SA with indices to bin
327: - index - selection of the index to bin on (SA_I or SA_J; PETSC_NULL is equivalent to SA_I)
330: Output Parameters:
331: . bins - SA containing the concatenated binned indices; the number of bins is the same as the result of ISGetImageSizeLocal().
333: Level: advanced
335: Concepts: binning^global indices
337: .seealso: SAMappingGetSupport(), SAMappingGetImage(), SAMappingGetSupportSizeLocal(), SAMappingGetImageSizeLocal(),
338: SAMappingBinLocal(), SAMappingMapLocal(), SAMappingMapLocal(), SAMappingBinSplit()
340: @*/
341: PetscErrorCode SAMappingBin(SAMapping map, SA array, SAIndex index, SA bins)
342: {
347: SAMappingCheckAssembled(map,PETSC_TRUE, 1);
348: if(!index) index = SA_I;
349: SAMappingCheckMethod(map,map->ops->bin,"SAMappingBin");
350: PetscLogEventBegin(SA_MAPPING_Bin,map,0,0,0);
351: (*map->ops->bin)(map,array,index,bins);
352: PetscLogEventEnd(SA_MAPPING_Bin,map,0,0,0);
353: return(0);
354: }
358: /*@C
359: SAMappingMapSplit - maps an SA with global indices from the rank's support to global indices from the rank's range.
360: The image of each index is a separate SA. See SAMappingMap, if concatenated output is desired.
361: Since SAMapping is in general multivalued, some global indices are mapped to multiple global indices.
362: Only the indices of the selected type (I or J) are mapped; the other indices and weights, if any,
363: are preserved on the images.
366: Not collective
368: Input Parameters:
369: + map - mapping of indices
370: . inarr - input SA
371: - index - selection of the index to map (SA_I or SA_J; PETSC_NULL is equivalent to SA_I)
374: Output Parameters:
375: . outarrs - SA list; the list length is the same as inarr's SA length.
378: Level: advanced
380: Concepts: mapping^indices global split
382: .seealso: SAMappingGetSupport(), SAMappingGetImage(), SAMappingGetSupportSizeLocal(), SAMappingGetImageSizeLocal(),
383: SAMappingMap(), SAMappingMapLocalSplit(), SAMappingBinSplit(), SAMappingBinSplitLocal()
385: @*/
386: PetscErrorCode SAMappingMapSplit(SAMapping map, SA inarr, SAIndex index, SA *outarr)
387: {
392: if(!index) index = SA_I;
393: SAMappingCheckAssembled(map,PETSC_TRUE,1);
394: SAMappingCheckMethod(map,map->ops->mapsplit,"SAMappingMapSplit");
395: PetscLogEventBegin(SA_MAPPING_MapSplit,map,0,0,0);
396: (*map->ops->mapsplit)(map,inarr,index,outarr);
397: PetscLogEventEnd(SA_MAPPING_MapSplit,map,0,0,0);
398: return(0);
399: }
403: /*@C
404: SAMappingMapSplitLocal - maps an SA with local indices from the rank's support to global indices from the rank's range.
405: The image of each index is a separate SA. Since SAMapping is in general multivalued, some local
406: indices are mapped to multiple global indices. Only the indices of the selected type (I or J) are mapped;
407: the other indices and weights, if any, are preserved on the images.
410: Not collective
412: Input Parameters:
413: + map - mapping of indices
414: . inarr - input SA
415: - index - selection of the index to map (SA_I or SA_J; PETSC_NULL is equivalent to SA_I)
418: Output Parameters:
419: . outarrs - SA list; the list length is the same as inarr's SA length.
422: Level: advanced
424: Concepts: mapping^indices local split
426: .seealso: SAMappingGetSupport(), SAMappingGetImage(), SAMappingGetSupportSizeLocal(), SAMappingGetImageSizeLocal(),
427: SAMappingMapLocal(), SAMappingMapSplit(), SAMappingBinSplit(), SAMappingBinSplitLocal()
429: @*/
430: PetscErrorCode SAMappingMapSplitLocal(SAMapping map, SA inarr, SAIndex index, SA *outarr)
431: {
436: if(!index) index = SA_I;
437: SAMappingCheckAssembled(map,PETSC_TRUE,1);
438: SAMappingCheckMethod(map,map->ops->mapsplitlocal,"SAMappingMapSplitLocal");
439: PetscLogEventBegin(SA_MAPPING_MapSplitLocal,map,0,0,0);
440: (*map->ops->mapsplitlocal)(map,inarr,index,outarr);
441: PetscLogEventEnd(SA_MAPPING_MapSplitLocal,map,0,0,0);
442: return(0);
443: }
447: /*@C
448: SAMappingBinSplitLocal - order local indices from the rank's support into n consecutive groups or "bins" (some possibly empty)
449: according to which of the n image indices they are mapped to on this rank. The bins are returned
450: as individual SAs. See SAMappingBinLocal() if concatenated bin output is desired.
451: Since SAMapping is potentially multivalued, the same index can appear in multiple bins.
452: The binning is done on the indices of the selected type (I or J); the other indices and weights, if any,
453: are moved to the appropriate bin together with the selected indices.
456: Not collective
458: Input Parameters:
459: + map - mapping of indices
460: . array - SA with indices to bin
461: - index - selection of the index to bin on (SA_I or SA_J; PETSC_NULL is equivalent to SA_I)
464: Output Parameters:
465: . bins - SA list of bins; the number of bins is the same as the result of ISGetImageSizeLocal().
467: Level: advanced
469: Concepts: binning^local indices split
471: .seealso: SAMappingGetSupport(), SAMappingGetImage(), SAMappingGetSupportSizeLocal(), SAMappingGetImageSizeLocal(),
472: SAMappingBinLocal(), SAMappingMapSplit(), SAMappingMapSplitLocal(), SAMappingBinSplit()
474: @*/
475: PetscErrorCode SAMappingBinSplitLocal(SAMapping map, SA array, SAIndex index, SA *bins)
476: {
481: SAMappingCheckAssembled(map,PETSC_TRUE, 1);
482: if(!index) index = SA_I;
483: SAMappingCheckMethod(map,map->ops->binsplitlocal,"SAMappingBinSplitLocal");
484: PetscLogEventBegin(SA_MAPPING_BinSplitLocal,map,0,0,0);
485: (*map->ops->binsplitlocal)(map,array,index,bins);
486: PetscLogEventEnd(SA_MAPPING_BinSplitLocal,map,0,0,0);
487: return(0);
488: }
493: /*@C
494: SAMappingBinSplit - group global indices from the rank's support into n groups or "bins" (some possibly empty)
495: according to which of the n image indices they are mapped to on this rank. The bins and
496: returned as individual SAs. See SAMappingBin() if concatenated bin output is desired.
497: Since SAMapping is potentially multivalued, the same index can appear in multiple bins.
498: The binning is done on the indices of selected type (I or J); the other indices and weights,
499: if any, are moved to the appropriate bin together with the selected indices.
502: Not collective
504: Input Parameters:
505: + map - mapping of indices
506: . array - SA with indices to bin
507: - index - selection of the index to bin on (SA_I or SA_J; PETSC_NULL is equivalent to SA_I)
510: Output Parameters:
511: . bins - SA list of bins; the number of bins is the same as the result of ISGetImageSizeLocal().
513: Level: advanced
515: Concepts: binning^global indices split
517: .seealso: SAMappingGetSupport(), SAMappingGetImage(), SAMappingGetSupportSizeLocal(), SAMappingGetImageSizeLocal(),
518: SAMappingBin(), SAMappingMapSplit(), SAMappingMapSplitLocal(), SAMappingBinSplitLocal()
520: @*/
521: PetscErrorCode SAMappingBinSplit(SAMapping map, SA array, SAIndex index, SA *bins)
522: {
527: SAMappingCheckAssembled(map,PETSC_TRUE, 1);
528: if(!index) index = SA_I;
529: SAMappingCheckMethod(map,map->ops->binsplit,"SAMappingBinSplit");
530: PetscLogEventBegin(SA_MAPPING_BinSplit,map,0,0,0);
531: (*map->ops->binsplit)(map,array,index,bins);
532: PetscLogEventEnd(SA_MAPPING_BinSplit,map,0,0,0);
533: return(0);
534: }
538: /*@C
539: SAMappingSetSizes - Sets the local and global sizes for the domain and range, and checks to determine compatibility
541: Collective on SAMapping
543: Input Parameters:
544: + map - the mapping
545: . m - number of local domain indices (or PETSC_DETERMINE)
546: . n - number of local range columns (or PETSC_DETERMINE)
547: . M - number of global domain indices (or PETSC_DETERMINE or PETSC_IGNORE)
548: - N - number of global range indices (or PETSC_DETERMINE or PETSC_IGNORE)
550: Notes:
551: The sizes specify (i) what the domain and range for the graph underlying SAMapping is,
552: and (b) what the parallel layout of the domain and range are. The local and global sizes
553: m and M (and n and N) are not independent. Two situations are possible the domain
554: (and,analogously, for the range):
555:
556: (P) Parallel layout for the domain demands that the local values of m on all ranks of
557: the communicator add up to M (see more at MatSetSizes, for example). Thus, m and M must
558: be specified in this compatible way. One way to ensure this is to specify m and leave
559: M as PETSC_DETERMINE -- then M is computed by summing the local m across the ranks.
560: The other option is to specify M (the same on all ranks, which will be checked) and
561: leave m as PETSC_DETERMINE. In this case the local m is determined by dividing M as
562: equally as possible among the ranks (m might end up being 0). If both m and M are specified,
563: their compatibility is verified by summing the m across the ranks. If m or M are PETSC_DETERMINE
564: on one rank, they must be PETSC_DETERMINE on all of the ranks, or the program might hang.
565: Finally, both m and M cannot be PETSC_DETERMINE at once.
567: In any case, domain indices can have any value 0 <= i < M on every rank (with the same M).
568: However, domain indices are split up among the ranks: each rank will "own" m indices with
569: the indices owned by rank 0 numbered [0,m), followed by the indices on rank 1, and so on.
571: (S) Sequential layout for the domain makes it essentially into a disjoint union of local
572: domains of local size m. This is signalled by specifying M as PETSC_IGNORE.
574: In this case, domain indices can have any value 0 <= i < m on every rank (with its own m).
576: Assembly/Mapping:
577: Whether the domain is laid out in parallel (P) or not (S), determines the behavior of SAMapping
578: during assembly. In case (P), the edges of the underlying graph are migrated to the rank that
579: owns the corresponding domain indices. SAMapping can map indices lying in its local range,
580: which is a subset of its local domain. This means that due to parallel assembly edges inserted
581: by different ranks might be used during the mapping. This is completely analogous to matrix
582: assembly.
584: When the domain is not laid out in parallel, no migration takes place and the mapping of indices
585: is done purely locally.
586:
587:
589: Support/Image:
590: Observe that the support and image of the graph may be strictly smaller than its domain and range,
591: if no edges from some domain points (or to some range points) are added to SAMapping.
592:
593: Operator:
594: Observe also that the linear operator defined by SAMapping will behave essentially as a VecScatter
595: (i) between MPI vectors with sizes (m,M) and (n,N), if both the domain and the range are (P),
596: (ii) between an MPI Vec with size (m,M) and a collection of SEQ Vecs (one per rank) of local size (n),
597: if the domain is (P) and the range is (S),
598: (iii) between a collection of SEQ Vecs (one per rank) of local size (m) and an MPI Vec of size (n,N),
599: if the domain is (S) and the range is (P),
600: (iv) between collections of SEQ Vecs (one per rank) or local sizes (m) and (n), if both the domain
601: and the range are (S).
603: Level: beginner
605: .seealso: SAMappingGetSizes(), SAMappingGetSupportSize(), SAMappingGetImageSize(), SAMappingMapIndicesLocal()
606: @*/
607: PetscErrorCode SAMappingSetSizes(SAMapping map, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
608: {
613: if (M > 0 && m > M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M);
614: if (N > 0 && n > N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N);
615: if(!map->xlayout) {
616: if(M == PETSC_IGNORE) {
617: PetscLayoutCreate(PETSC_COMM_SELF, &(map->xlayout));
618: }
619: else {
620: PetscLayoutCreate(((PetscObject)map)->comm, &(map->xlayout));
621: }
622: }
623: if(!map->ylayout) {
624: if(N == PETSC_IGNORE) {
625: PetscLayoutCreate(PETSC_COMM_SELF, &(map->ylayout));
626: }
627: else {
628: PetscLayoutCreate(((PetscObject)map)->comm, &(map->ylayout));
629: }
630: }
631: if ((map->xlayout->n >= 0 || map->xlayout->N >= 0) && (map->xlayout->n != m || map->xlayout->N != M)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset domain sizes to %D local %D global after previously setting them to %D local %D global",m,M,map->xlayout->n,map->xlayout->N);
632: if ((map->ylayout->n >= 0 || map->ylayout->N >= 0) && (map->ylayout->n != n || map->ylayout->N != N)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset range sizes to %D local %D global after previously setting them to %D local %D global",n,N,map->ylayout->n,map->ylayout->N);
633:
634: map->xlayout->n = m;
635: map->ylayout->n = n;
636: map->xlayout->N = M;
637: map->ylayout->N = N;
639: map->setup = PETSC_FALSE;
640: map->assembled = PETSC_FALSE;
641: return(0);
642: }
646: PetscErrorCode SAMappingGetSizes(SAMapping map, PetscInt *m, PetscInt *n, PetscInt *M, PetscInt *N)
647: {
650: *m = map->xlayout->n;
651: *n = map->ylayout->n;
652: *M = map->xlayout->N;
653: *N = map->ylayout->N;
654: return(0);
655: }
659: PetscErrorCode SAMappingSetUp_SAMapping(SAMapping map)
660: {
663: PetscLayoutSetBlockSize(map->xlayout,1);
664: PetscLayoutSetBlockSize(map->ylayout,1);
665: PetscLayoutSetUp(map->xlayout);
666: PetscLayoutSetUp(map->ylayout);
667: return(0);
668: }/* SAMappingSetUp_SAMapping() */
672: /*@C
673: SAMappingSetUp - Sets up the internal mapping data structures for the later use.
675: Collective on SAMapping
677: Input Parameters:
678: . map - the SAMapping context
680: Notes:
681: For basic use of the SAMapping classes the user need not explicitly call
682: SAMappingSetUp(), since these actions will happen automatically.
684: Level: advanced
686: .keywords: SAMapping, setup
688: .seealso: SAMappingCreate(), SAMappingDestroy(), SAMappingSetSizes()
689: @*/
690: PetscErrorCode SAMappingSetUp(SAMapping map)
691: {
695: SAMappingCheckMethod(map,map->ops->setup,"SAMappingSetUp");
696: (*(map->ops->setup))(map);
697: map->setup = PETSC_TRUE;
698: return(0);
699: }/* SAMappingGetSizes() */
701: #undef __FUNCT__
703: PetscErrorCode SAMappingAssemblyBegin(SAMapping map)
704: {
709: if(map->assembled) return(0);
710: SAMappingSetUp(map);
711:
712: SAMappingCheckMethod(map,map->ops->assemblybegin, "SAMappingAsemblyBegin");
713: PetscLogEventBegin(SA_MAPPING_AssemblyBegin, map,0,0,0);
714: (*(map->ops->assemblybegin))(map);
715: PetscLogEventEnd(SA_MAPPING_AssemblyBegin, map,0,0,0);
717: return(0);
718: }/* SAMappingAssemblyBegin() */
722: PetscErrorCode SAMappingAssemblyEnd(SAMapping map)
723: {
727: SAMappingCheckMethod(map,map->ops->assemblyend, "SAMappingAsemblyEnd");
728: PetscLogEventBegin(SA_MAPPING_AssemblyEnd, map,0,0,0);
729: (*(map->ops->assemblyend))(map);
730: PetscLogEventEnd(SA_MAPPING_AssemblyBegin, map,0,0,0);
731: map->assembled = PETSC_TRUE;
732: return(0);
733: }/* SAMappingAssemblyEnd() */
736: #undef __FUNCT__
738: PetscErrorCode SAMappingGetSupport(SAMapping map, PetscInt *_len, PetscInt *_supp[]) {
741:
743: SAMappingCheckAssembled(map,PETSC_TRUE,1);
744: if(!_len && !_supp) return(0);
745: SAMappingCheckMethod(map,map->ops->getsupport,"SAMappingGetSupport");
746: (*(map->ops->getsupport))(map,_len,_supp);
747: return(0);
748: }
750: #undef __FUNCT__
752: PetscErrorCode SAMappingGetSupportIS(SAMapping map, IS *supp) {
755:
757: SAMappingCheckAssembled(map,PETSC_TRUE,1);
758: if(!supp) return(0);
759: SAMappingCheckMethod(map,map->ops->getsupportis,"SAMappingGetSupportIS");
760: (*(map->ops->getsupportis))(map,supp);
761: return(0);
762: }
765: #undef __FUNCT__
767: PetscErrorCode SAMappingGetSupportSA(SAMapping map, SA *supp) {
770:
772: SAMappingCheckAssembled(map,PETSC_TRUE,1);
773: if(!supp) return(0);
774: SAMappingCheckMethod(map,map->ops->getsupportsa,"SAMappingGetSupportSA");
775: (*(map->ops->getsupportsa))(map,supp);
776: return(0);
777: }
779: #undef __FUNCT__
781: PetscErrorCode SAMappingGetImage(SAMapping map, PetscInt *_len, PetscInt *_image[]) {
784:
786: SAMappingCheckAssembled(map,PETSC_TRUE,1);
787: if(!_len && !_image) return(0);
788: SAMappingCheckMethod(map,map->ops->getimage,"SAMappingGetImage");
789: (*(map->ops->getimage))(map,_len,_image);
790: return(0);
791: }
793: #undef __FUNCT__
795: PetscErrorCode SAMappingGetImageIS(SAMapping map, IS *image) {
798:
800: SAMappingCheckAssembled(map,PETSC_TRUE,1);
801: if(!image) return(0);
802: SAMappingCheckMethod(map,map->ops->getimageis,"SAMappingGetImageIS");
803: (*(map->ops->getimageis))(map,image);
804: return(0);
805: }
807: #undef __FUNCT__
809: PetscErrorCode SAMappingGetImageSA(SAMapping map, SA *image) {
812:
814: SAMappingCheckAssembled(map,PETSC_TRUE,1);
815: if(!image) return(0);
816: SAMappingCheckMethod(map,map->ops->getimagesa,"SAMappingGetImageSA");
817: (*(map->ops->getimagesa))(map,image);
818: return(0);
819: }
823: PetscErrorCode SAMappingGetMaxImageSize(SAMapping map, PetscInt *maxsize)
824: {
829: SAMappingCheckAssembled(map,PETSC_TRUE,1);
830: SAMappingCheckMethod(map, map->ops->getmaximagesize,"SAMappingGetMaxImageSize");
831: (*map->ops->getmaximagesize)(map,maxsize);
832: return(0);
833: }
838: PetscErrorCode SAMappingGetOperator(SAMapping map, Mat *mat)
839: {
844: SAMappingCheckAssembled(map,PETSC_TRUE,1);
845: SAMappingCheckMethod(map, map->ops->getoperator,"SAMappingGetOperator");
846: (*map->ops->getoperator)(map,mat);
848: return(0);
849: }
851: #undef __FUNCT__
853: PetscErrorCode SAMappingView(SAMapping map, PetscViewer v)
854: {
858: SAMappingCheckAssembled(map,PETSC_TRUE,1);
859: SAMappingCheckMethod(map,map->ops->view, "SAMappingView");
860: (*(map->ops->view))(map,v);
861: return(0);
862: }/* SAMappingView() */
864: #undef __FUNCT__
866: PetscErrorCode SAMappingInvert(SAMapping map, SAMapping *imap)
867: {
872: SAMappingCheckMethod(map,map->ops->invert, "SAMappingInvert");
873: PetscLogEventBegin(SA_MAPPING_Invert, map, 0,0,0);
874: (*(map->ops->invert))(map,imap);
875: PetscLogEventEnd(SA_MAPPING_Invert, map, 0,0,0);
876: return(0);
877: }
881: /*@C
882: SAMappingPullback - compose mappings C = A*B
884: Collective on SAMapping
886: Input Parameters:
887: + A - the left mapping
888: - B - the right mapping
891: Output Parameters:
892: . C - the product mapping: domain as in A, range as in B
895: Level: intermediate
897: .seealso: SAMappingPushforward(), SAMappingInvert()
898: @*/
899: PetscErrorCode SAMappingPullback(SAMapping A,SAMapping B, SAMapping *C)
900: {
902: PetscErrorCode (*pullback)(SAMapping,SAMapping,SAMapping*);
903: char pullbackname[256];
908: if (!A->assembled) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled mapping");
911: if (!B->assembled) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled mapping");
914: /* dispatch based on the type of A and B */
915: PetscStrcpy(pullbackname,"SAMappingPullback_");
916: PetscStrcat(pullbackname,((PetscObject)A)->type_name);
917: PetscStrcat(pullbackname,"_");
918: PetscStrcat(pullbackname,((PetscObject)B)->type_name);
919: PetscStrcat(pullbackname,"_C"); /* e.g., pullbackname = "ISPullback_ismappingis_ismappingis_C" */
920: PetscObjectQueryFunction((PetscObject)B,pullbackname,(void (**)(void))&pullback);
921: if (!pullback) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"SAMappingPullback requires A, %s, to be compatible with B, %s",((PetscObject)A)->type_name,((PetscObject)B)->type_name);
922: PetscLogEventBegin(SA_MAPPING_Pullback, A,B,0,0);
923: (*pullback)(A,B,C);
924: PetscLogEventEnd(SA_MAPPING_Pullback, A,B,0,0);
925: return(0);
926: }
930: /*@C
931: SAMappingPushforward - mapping from the range of A to the range of B, pointwise on the common domain
933: Collective on SAMapping
935: Input Parameters:
936: + A - the left mapping
937: - B - the right mapping
940: Output Parameters:
941: . C - the product mapping: domain as the range of A, range as the range of B
944: Level: intermediate
945: .seealso: SAMappingPullback(), SAMappingInvert()
946: @*/
947: PetscErrorCode SAMappingPushforward(SAMapping A,SAMapping B, SAMapping *C)
948: {
950: PetscErrorCode (*pushforward)(SAMapping,SAMapping,SAMapping*);
951: char pushforwardname[256];
956: if (!A->assembled) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled mapping");
959: if (!B->assembled) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled mapping");
962: /* dispatch based on the type of A and B */
963: PetscStrcpy(pushforwardname,"SAMappingPushforward_");
964: PetscStrcat(pushforwardname,((PetscObject)A)->type_name);
965: PetscStrcat(pushforwardname,"_");
966: PetscStrcat(pushforwardname,((PetscObject)B)->type_name);
967: PetscStrcat(pushforwardname,"_C"); /* e.g., pushforwardname = "ISPushforward_ismappingis_ismappingis_C" */
968: PetscObjectQueryFunction((PetscObject)B,pushforwardname,(void (**)(void))&pushforward);
969: if (!pushforward) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"SAMappingPushforward requires A, %s, to be compatible with B, %s",((PetscObject)A)->type_name,((PetscObject)B)->type_name);
970: PetscLogEventBegin(SA_MAPPING_Pushforward, A,B,0,0);
971: (*pushforward)(A,B,C);
972: PetscLogEventEnd(SA_MAPPING_Pushforward, A,B,0,0);
973: return(0);
974: }
976: #undef __FUNCT__
978: PetscErrorCode SAMappingSetType(SAMapping map, const SAMappingType maptype) {
980: PetscErrorCode (*ctor)(SAMapping);
981: PetscBool sametype;
984: PetscTypeCompare((PetscObject)map,maptype, &sametype);
985: if(sametype) return(0);
987: if(!SAMappingRegisterAllCalled) {
988: SAMappingRegisterAll(PETSC_NULL);
989: }
990: PetscFListFind(ISList,((PetscObject)map)->comm,maptype,PETSC_TRUE,(void(**)(void))&ctor);
991: if(!ctor) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized SAMapping type: %s", maptype);
993: /* destroy the old implementation, if it existed */
994: if(map->ops->destroy) {
995: (*(map->ops->destroy))(map);
996: map->ops->destroy = PETSC_NULL;
997: }
998:
999: /* create the new implementation */
1000: (*ctor)(map);
1001: return(0);
1002: }/* SAMappingSetType() */
1004: #undef __FUNCT__
1006: PetscErrorCode SAMappingDestroy(SAMapping *_map)
1007: {
1010: if(!*_map) return(0);
1012: if(--((PetscObject)(*_map))->refct > 0) return(0);
1013: if((*_map)->ops->destroy) {
1014: (*(*_map)->ops->destroy)(*_map);
1015: }
1016: PetscLayoutDestroy(&(*_map)->xlayout);
1017: PetscLayoutDestroy(&(*_map)->ylayout);
1018: PetscHeaderDestroy(_map);
1019: return(0);
1020: }/* SAMappingDestroy() */
1022: #undef __FUNCT__
1024: PetscErrorCode SAMappingCreate(MPI_Comm comm, SAMapping *_map)
1025: {
1026: SAMapping map;
1030: *_map = PETSC_NULL;
1031: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
1032: SAMappingInitializePackage(PETSC_NULL);
1033: #endif
1035: PetscHeaderCreate(map,_p_SAMapping,struct _SAMappingOps,SA_MAPPING_CLASSID,0,"SAMapping","SAMapping","DM",comm,SAMappingDestroy,SAMappingView);
1036: PetscLayoutCreate(comm,&(map->xlayout));
1037: PetscLayoutCreate(comm,&(map->ylayout));
1038: *_map = map;
1039: return(0);
1040: }/* SAMappingCreate() */
1044: #undef __FUNCT__
1046: PetscErrorCode SAHunkCreate(PetscInt maxlength, SAComponents mask, SAHunk *_hunk)
1047: {
1049: SAHunk hunk;
1053: if(maxlength <= 0)
1054: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Nonpositive SAHunk maxlength: %D", maxlength);
1056: PetscNew(struct _n_SAHunk, &hunk);
1057: hunk->mask = mask;
1058: hunk->maxlength = maxlength;
1059: if(mask & SA_I) {
1060: PetscInt *ia;
1061: PetscMalloc(sizeof(PetscInt)*maxlength, &ia);
1062: PetscMemzero(ia, sizeof(PetscInt)*maxlength);
1063: hunk->i = ia;
1064: }
1065: if(mask & SA_J) {
1066: PetscInt *ja;
1067: PetscMalloc(sizeof(PetscInt)*maxlength, &ja);
1068: PetscMemzero(ja, sizeof(PetscInt)*maxlength);
1069: hunk->j = ja;
1070: }
1071: if(mask & SA_W) {
1072: PetscScalar *wa;
1073: PetscMalloc(sizeof(PetscScalar)*maxlength, &wa);
1074: PetscMemzero(wa, sizeof(PetscScalar)*maxlength);
1075: hunk->w = wa;
1076: }
1077: hunk->mode = PETSC_OWN_POINTER;
1078: hunk->length = 0;
1079: hunk->refcnt = 1;
1080: *_hunk = hunk;
1081: return(0);
1082: }
1084: #undef __FUNCT__
1086: PetscErrorCode SAHunkAddData(SAHunk hunk, PetscInt length, const PetscInt *i, const PetscScalar *w, const PetscInt *j)
1087: {
1089: PetscInt mask;
1092: mask = (i != PETSC_NULL) | ((j != PETSC_NULL)<<1) | ((w != PETSC_NULL)<<2);
1093: if(mask != hunk->mask) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data components %D incompatible with the SAHunk mask", mask,hunk->mask);
1094: if(hunk->length + length > hunk->maxlength) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot add data of length %D, hunk only has %D space left", length, hunk->maxlength-hunk->length);
1095: if(mask & SA_I) {
1096: PetscMemcpy(hunk->i+hunk->length, i, sizeof(PetscInt)*length);
1097: }
1098: if(mask & SA_J) {
1099: PetscMemcpy(hunk->j+hunk->length, j, sizeof(PetscInt)*length);
1100: }
1101: if(mask & SA_W) {
1102: PetscMemcpy(hunk->w+hunk->length, w, sizeof(PetscScalar)*length);
1103: }
1104: hunk->length += length;
1105: return(0);
1106: }
1109: #undef __FUNCT__
1111: PetscErrorCode SAHunkGetSubHunk(SAHunk hunk, PetscInt start, PetscInt maxlength, PetscInt length, SAComponents mask, SAHunk *_subhunk)
1112: {
1114: SAHunk subhunk;
1118: *_subhunk = PETSC_NULL;
1119: if(maxlength <= 0)
1120: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Nonpositive subhunk maxlength: %D", maxlength);
1121: if(length <= 0)
1122: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Nonpositive subhunk length: %D", length);
1123: if(length > maxlength)
1124: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "subhunk length %D exceeds maxlength %D", length, maxlength);
1125: if(start < 0 || start >= hunk->maxlength)
1126: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "subhunk start %D out of bounds [0,%D)", start,hunk->maxlength);
1127: if(start+maxlength > hunk->maxlength)
1128: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "subhunk end %D beyond the end of hunk %D",start+maxlength,hunk->maxlength);
1129:
1130: if(mask & (~(hunk->mask)))
1131: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subhunk mask %D is not a submask of the hunk mask %D", mask, hunk->mask);
1132: PetscNew(struct _n_SAHunk, &subhunk);
1133: subhunk->mask = mask;
1134: subhunk->maxlength = maxlength;
1135: subhunk->length = length;
1136: if(mask & SA_I) {
1137: subhunk->i = hunk->i+start;
1138: }
1139: if(mask & SA_J) {
1140: subhunk->j = hunk->j+start;
1141: }
1142: if(mask & SA_W) {
1143: subhunk->w = hunk->w+start;
1144: }
1145: subhunk->mode = PETSC_USE_POINTER;
1146: ++(hunk->refcnt);
1147: *_subhunk = subhunk;
1148: return(0);
1149: }
1152: #undef __FUNCT__
1154: PetscErrorCode SAHunkDestroy(SAHunk *_hunk)
1155: {
1158: if(!*_hunk) return(0);
1159: if(--((*_hunk)->refcnt) > 0) {*_hunk = 0; return(0);}
1160: if((*_hunk)->length && (*_hunk)->mode != PETSC_USE_POINTER) {
1161: PetscFree((*_hunk)->i);
1162: PetscFree((*_hunk)->w);
1163: PetscFree((*_hunk)->j);
1164: }
1165: (*_hunk)->length = 0;
1166: if((*_hunk)->parent) {
1167: SAHunkDestroy(&((*_hunk)->parent));
1168: }
1169: PetscFree((*_hunk));
1170: *_hunk = PETSC_NULL;
1171: return(0);
1172: }
1175: #undef __FUNCT__
1177: PetscErrorCode SACreate(SAComponents mask, SA *_arr)
1178: {
1182: if(!(mask & SA_I) && !(mask & SA_J)) {
1183: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "SAComponents %D must contain at least one of the indices: I or J", mask);
1184: }
1185: PetscNew(struct _n_SA, _arr);
1186: (*_arr)->mask = mask;
1187: return(0);
1188: }
1190: #undef __FUNCT__
1192: PetscErrorCode SADuplicate(SA arr, SA *_darr)
1193: {
1195: SA darr;
1196: SALink link;
1200: SACreate(arr->mask, &darr);
1201: link = arr->first;
1202: while(link) {
1203: SAAddHunk(darr,link->hunk);
1204: }
1205: *_darr = arr;
1206: return(0);
1207: }
1209: #undef __FUNCT__
1211: PetscErrorCode SACreateArrays(SAComponents mask, PetscInt count, SA **_arrays)
1212: {
1214: PetscInt i;
1217: if(!(mask & SA_I) && !(mask & SA_J)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "SAComponents %D must contain at least one of the indices: I or J", mask);
1218: if(count < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Negative array count: %D", count);
1219: PetscMalloc(sizeof(SA), _arrays);
1220: for(i = 0; i < count; ++i) {
1221: SACreate(mask, *_arrays+i);
1222: }
1223: return(0);
1224: }
1226: #undef __FUNCT__
1228: PetscErrorCode SAClear(SA arr)
1229: {
1231: SALink link, next;
1234: link = arr->first;
1235: while(link) {
1236: SAHunkDestroy(&(link->hunk));
1237: next = link->next;
1238: PetscFree(link);
1239: link = next;
1240: }
1241: arr->first = PETSC_NULL;
1242: arr->last = PETSC_NULL;
1243: arr->length = 0;
1244: return(0);
1245: }
1248: #undef __FUNCT__
1250: PetscErrorCode SADestroy(SA *_arr)
1251: {
1254: if(!*_arr) return(0);
1255: SAClear(*_arr);
1256: (*_arr)->mask = 0;
1257: PetscFree(*_arr);
1258: *_arr = PETSC_NULL;
1259: return(0);
1260: }
1263: /*
1264: Right now we merely allocate a new hunk.
1265: In the future, this can manage a pool of hunks, use a buffer to draw subhunks from, etc.
1266: */
1267: #undef __FUNCT__
1269: PetscErrorCode SAGetHunk(SA chain, PetscInt length, SAHunk *_hunk)
1270: {
1273: SAHunkCreate(length, chain->mask, _hunk);
1274: return(0);
1275: }
1277: #undef __FUNCT__
1279: PetscErrorCode SAAddHunk(SA chain, SAHunk hunk)
1280: {
1282: SALink link;
1284: if(chain->mask & (~(hunk->mask))) {
1285: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Hunk mask %D incompatible with the array mask %D", hunk->mask, chain->mask);
1286: }
1287: PetscMalloc(sizeof(struct _n_SALink), &link);
1288: link->hunk = hunk;
1289: ++(hunk->refcnt);
1290: if(chain->last) {
1291: chain->last->next = link;
1292: chain->last = link;
1293: }
1294: else {
1295: chain->first = chain->last = link;
1296: }
1297: return(0);
1298: }
1301: #undef __FUNCT__
1303: PetscErrorCode SAAddData(SA chain, PetscInt length, const PetscInt *i, const PetscScalar *w, const PetscInt *j)
1304: {
1306: SAHunk hunk;
1307: SAComponents mask;
1310: mask = (i != PETSC_NULL) | ((j != PETSC_NULL)<<1) | ((w != PETSC_NULL)<<2);
1311: if(mask != chain->mask) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data components %D incompatible with the SAComponents", mask,chain->mask);
1312: SAGetHunk(chain, length, &hunk);
1313: SAHunkAddData(hunk, length, i,w,j);
1314: SAAddHunk(chain, hunk);
1315: return(0);
1316: }
1319: #undef __FUNCT__
1321: PetscErrorCode SAAddI(SA chain, PetscInt length, PetscInt i, const PetscScalar wa[], const PetscInt ja[])
1322: {
1324: SAHunk hunk;
1325: PetscInt mask;
1326: PetscInt k;
1329: mask = (SA_I | (ja != PETSC_NULL)<<1 | (wa != PETSC_NULL)<<2);
1330: if(mask & (~(chain->mask))) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Array components provided %D incompatible with array mask %D", mask, chain->mask);
1331: if(!length) return(0);
1332: SAGetHunk(chain, length, &hunk);
1333: for(k = 0; k < length; ++k) hunk->i[k] = i;
1334: if(ja) {
1335: PetscMemcpy(hunk->j, ja, sizeof(PetscInt)*length);
1336: }
1337: if(wa) {
1338: PetscMemcpy(hunk->w, wa, sizeof(PetscScalar)*length);
1339: }
1340: SAAddHunk(chain, hunk);
1341: return(0);
1342: }
1344: #undef __FUNCT__
1346: PetscErrorCode SAAddJ(SA chain, PetscInt length, const PetscInt ia[], const PetscScalar wa[], PetscInt j)
1347: {
1349: SAHunk hunk;
1350: PetscInt mask;
1351: PetscInt k;
1354: mask = ((PetscInt)(ia != PETSC_NULL) | (PetscInt)SA_J | (PetscInt)(wa != PETSC_NULL)<<2);
1355: if(mask & (~(chain->mask))) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Array components provided %D incompatible with array mask %D", mask, chain->mask);
1356: if(!length) return(0);
1357: SAGetHunk(chain, length, &hunk);
1358: for(k = 0; k < length; ++k) hunk->j[k] = j;
1359: if(ia) {
1360: PetscMemcpy(hunk->i, ia, sizeof(PetscInt)*length);
1361: }
1362: if(wa) {
1363: PetscMemcpy(hunk->w, wa, sizeof(PetscScalar)*length);
1364: }
1365: SAAddHunk(chain, hunk);
1366: return(0);
1367: }
1369: #undef __FUNCT__
1371: PetscErrorCode SAMerge(SA arr, SA *_marr)
1372: {
1374: SALink link;
1375: SAHunk merged;
1376: PetscInt count, offset;
1381: /* Determine the number of links in the chain and perform the mask consistency check. */
1382: link = arr->first;
1383: count = 0;
1384: while(link) {
1385: /* Mask consistency check. */
1386: if(arr->mask & (~(link->hunk->mask))) {
1387: SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Hunk components %D in link %D incompatible with the merge mask %D", link->hunk->mask, count, arr->mask);
1388: }
1389: link = link->next;
1390: ++count;
1391: }
1392: if(count == 1) return(0);
1394: if(arr->length) {
1395: SAHunkCreate(arr->length, arr->mask, &merged);
1396: /* Copy the indices and weights into the merged arrays. */
1397: offset = 0;
1398: link = arr->first;
1399: while(link) {
1400: SAHunk hunk = link->hunk;
1401: if(arr->mask & SA_I) {
1402: PetscMemcpy(merged->i+offset, hunk->i, sizeof(PetscInt)*hunk->length);
1403: }
1404: if(arr->mask & SA_J) {
1405: PetscMemcpy(merged->j+offset, hunk->j, sizeof(PetscInt)*hunk->length);
1406: }
1407: if(arr->mask & SA_W) {
1408: PetscMemcpy(merged->w+offset, hunk->w, sizeof(PetscScalar)*hunk->length);
1409: }
1410: offset += hunk->length;
1411: }
1412: }/* if(arr->length) */
1413: merged->length = offset;
1414: SACreate(arr->mask, _marr);
1415: SAAddHunk(*_marr, merged);
1416: return(0);
1417: }
1419: #undef __FUNCT__
1421: PetscErrorCode SAGetLength(SA chain, PetscInt *_length)
1422: {
1426: *_length = chain->length;
1427: return(0);
1428: }
1430: #undef __FUNCT__
1432: PetscErrorCode SAGetData(SA chain, PetscInt *ia, PetscScalar *wa, PetscInt *ja)
1433: {
1435: PetscInt mask, off;
1436: SALink link;
1437: SAHunk hunk;
1441: mask = ((ia != PETSC_NULL) & (chain->mask&SA_I)) | ((ja != PETSC_NULL)<<1 & (chain->mask&SA_J)) | ((wa != PETSC_NULL)<<2 & (chain->mask&SA_W));
1442: off = 0;
1443: link = chain->first;
1444: while(link) {
1445: hunk = link->hunk;
1446: if(mask&SA_I) {
1447: PetscMemcpy(ia+off, hunk->i, sizeof(PetscInt)*hunk->length);
1448: }
1449: if(mask&SA_J) {
1450: PetscMemcpy(ja+off, hunk->j, sizeof(PetscInt)*hunk->length);
1451: }
1452: if(mask&SA_W) {
1453: PetscMemcpy(wa+off, hunk->w, sizeof(PetscScalar)*hunk->length);
1454: }
1455: link = link->next;
1456: off += hunk->length;
1457: }
1458: return(0);
1459: }
1461: #undef __FUNCT__
1463: PetscErrorCode SAAddArray(SA chain, SA chain2)
1464: {
1466: SALink link;
1470: if(chain->mask & (~(chain2->mask))) {
1471: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "New mask %D ncompatible with original %D", chain2->mask, chain->mask);
1472: }
1473: link = chain2->first;
1474: while(link) {
1475: SAAddHunk(chain, link->hunk);
1476: link = link->next;
1477: }
1478: return(0);
1479: }
1481: #undef __FUNCT__
1483: PetscErrorCode SASplit(SA arr, PetscInt count, const PetscInt *lengths, PetscInt mask, SA *arrs)
1484: {
1486: PetscInt i, lengthi, start, end, len;
1487: SALink link;
1488: SAHunk hunk, subhunk;
1492: if(count < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Negative array length %D", count);
1495: if(!count) return(0);
1496: link = arr->first;
1497: i = 0; /* current subarray index */
1498: lengthi = 0; /* current length of the current subarray */
1499: while(link) {
1500: hunk = link->hunk;
1501: start = 0; /* start of unassigned hunk space */
1502: while(hunk->length-start) {
1503: end = PetscMin(hunk->length, lengths[i]-lengthi);
1504: len = start-end;
1505: SAHunkGetSubHunk(hunk,start,len,len,hunk->mask, &subhunk);
1506: SAAddHunk(arrs[i], subhunk);
1507: start += len;
1508: lengthi += len;
1509: if(lengthi == lengths[i]) {
1510: lengthi = 0; ++i;
1511: if(i >= count) return(0);
1512: }
1513: }
1514: link = link->next;
1515: }
1516:
1517: return(0);
1518: }
1522: /*
1523: Checks whether all indices are within [imin,imax) and generate an error, if they are not and
1524: if outOfBoundsError == PETSC_TRUE. Return the result in flag.
1525: */
1529: {
1530: PetscInt i;
1531: PetscBool inBounds = PETSC_TRUE;
1533: for (i=0; i<len; ++i) {
1534: if (idx[i] < imin) {
1535: if(outOfBoundsError) {
1536: SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D at location %D is less than min %D",idx[i],i,imin);
1537: }
1538: inBounds = PETSC_FALSE;
1539: break;
1540: }
1541: if (idx[i] >= imax) {
1542: if(outOfBoundsError) {
1543: SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D at location %D is greater than max %D",idx[i],i,imax);
1544: }
1545: inBounds = PETSC_FALSE;
1546: break;
1547: }
1548: }
1549: if(flag) *flag = inBounds;
1550: return(0);
1551: }
1553: /*
1554: Checks if any indices are within [imin,imax) and generate an error, if they are not and
1555: if outOfBoundsError == PETSC_TRUE. Return the result in flag.
1556: */
1560: {
1561: PetscInt n;
1562: PetscBool inBounds = PETSC_TRUE, isstride;
1565:
1566: ISGetLocalSize(is, &n);
1567: PetscTypeCompare((PetscObject)is,ISSTRIDE,&isstride);
1568: if(isstride) {
1569: PetscInt first, step, last;
1570:
1571: ISStrideGetInfo(is, &first, &step);
1572: last = first + step*n;
1573: if (first < imin || last < imin) {
1574: if(outOfBoundsError)
1575: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index smaller than min %D", imin);
1576: inBounds = PETSC_FALSE;
1577: goto functionend;
1578: }
1579: if (first >= imax || last >= imax) {
1580: if(outOfBoundsError)
1581: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index greater than max %D", imax);
1582: inBounds = PETSC_FALSE;
1583: goto functionend;
1584: }
1585: } else { /* not stride */
1586: const PetscInt *idx;
1587: ISGetIndices(is, &idx);
1589: ISRestoreIndices(is, &idx);
1590: }
1591: functionend:
1592: if(flag) *flag = inBounds;
1593: return(0);
1594: }
1596: /*
1597: FIX: If we wanted to split this into AssemblyBegin/End, some data must be passed between the stages (e.g., tags,
1598: waits), this could be made into an SAMapping method. However, then an SAMapping of some sort must be instantiated
1599: for SA assembly. Maybe it always happens in the common use cases, anyhow.
1600: */
1603: PetscErrorCode SAAssemble(SA chain, PetscInt mask, PetscLayout layout, SA achain)
1604: {
1606: MPI_Comm comm = layout->comm;
1607: PetscMPIInt size, rank, tag_i, tag_v;
1608: PetscMPIInt chainmaskmpi, allchainmask;
1609: const PetscInt *ixidx, *iyidx;
1610: PetscInt idx, lastidx;
1611: PetscInt i, j, p;
1612: PetscInt *owner = PETSC_NULL;
1613: PetscMPIInt nsends, nrecvs;
1614: PetscMPIInt *plengths, *sstarts = PETSC_NULL;
1615: PetscMPIInt *rnodes, *rlengths, *rstarts = PETSC_NULL, rlengthtotal;
1616: PetscInt **rindices = PETSC_NULL, *sindices= PETSC_NULL;
1617: PetscScalar *svalues = PETSC_NULL, **rvalues = PETSC_NULL;
1618: MPI_Request *recv_reqs_i, *recv_reqs_v, *send_reqs;
1619: MPI_Status recv_status, *send_statuses;
1620: #if defined(PETSC_USE_DEBUG)
1621: PetscBool found;
1622: #endif
1623: PetscInt ni, nv, count, alength;
1624: PetscInt *aixidx, *aiyidx = PETSC_NULL;
1625: PetscScalar *aval = PETSC_NULL;
1626: const PetscScalar *val;
1627: SALink link;
1628: SAHunk hunk, ahunk;
1634: /* Make sure that at least one of the indices -- I or J -- is being assembled on, and that the index being assembled on is present in the SA. */
1635: if((mask != SA_I && mask != SA_J)|| !(mask & chain->mask)) {
1636: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot assemble SA with components %D on component %D", chain->mask, mask);
1637: }
1638: if(chain->mask != achain->mask) {
1639: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input and output array masks differ: %D and %D", chain->mask, achain->mask);
1640: }
1642: /* Comm parameters */
1643: MPI_Comm_size(comm,&size);
1644: MPI_Comm_rank(comm,&rank);
1646: /* If layout isn't parallel, then this is a noop. */
1647: if(size == 1) return(0);
1649: /* Make sure that chain type is the same across the comm. */
1650: chainmaskmpi = PetscMPIIntCast(chain->mask);
1651: MPI_Allreduce(&(chainmaskmpi), &allchainmask, 1, MPI_INT, MPI_BOR, comm);
1652: if(allchainmask^chainmaskmpi) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Chain mask must be the same across the communicator. Got %D globally and %D locally", allchainmask, chainmaskmpi);
1654: /* How many index arrays are being sent? One or two? */
1655: ni = 0;
1656: ni += ((chain->mask & SA_I) > 0);
1657: ni += ((chain->mask & SA_J) > 0);
1659: /* How many value arrays are being sent? One or none? */
1660: nv = ((chain->mask & SA_W) > 0);
1663:
1664: /*
1665: Each processor ships off its ixidx[j] and, possibly, the appropriate combination of iyidx[j]
1666: and val[j] to the appropriate processor.
1667: */
1668: /* first count number of contributors to each processor */
1669: PetscMalloc2(size,PetscMPIInt,&plengths,chain->length,PetscInt,&owner);
1670: PetscMemzero(plengths,size*sizeof(PetscMPIInt));
1671: lastidx = -1;
1672: count = 0;
1673: p = 0;
1674: link = chain->first;
1675: while(link) {
1676: hunk = link->hunk;
1677: for (i=0; i<hunk->length; ++i) {
1678: if(mask == SA_I) {
1679: ixidx = hunk->i;
1680: iyidx = hunk->j;
1681: }
1682: else {
1683: ixidx = hunk->j;
1684: iyidx = hunk->i;
1685: }
1686: /* if indices are NOT locally sorted, need to start search for the proc owning inidx[i] at the beginning */
1687: if (lastidx > (idx = ixidx[i])) p = 0;
1688: lastidx = idx;
1689: for (; p<size; ++p) {
1690: if (idx >= layout->range[p] && idx < layout->range[p+1]) {
1691: plengths[p]++;
1692: owner[count] = p;
1693: #if defined(PETSC_USE_DEBUG)
1694: found = PETSC_TRUE;
1695: #endif
1696: break;
1697: }
1698: }
1699: #if defined(PETSC_USE_DEBUG)
1700: if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D out of range",idx);
1701: found = PETSC_FALSE;
1702: #endif
1703: ++count;
1704: }/* for(i=0; i < hunk->length; ++i) */
1705: link = link->next;
1706: }
1707: nsends = 0; for (p=0; p<size; ++p) { nsends += (plengths[p] > 0);}
1708:
1709: /* inform other processors of number of messages and max length*/
1710: PetscGatherNumberOfMessages(comm,PETSC_NULL,plengths,&nrecvs);
1711: PetscGatherMessageLengths(comm,nsends,nrecvs,plengths,&rnodes,&rlengths);
1712: /* no values are sent if nv == 0 */
1713: if(nv) {
1714: /* values in message i are rlengths[i] in number */
1715: PetscCommGetNewTag(layout->comm, &tag_v);
1716: PetscPostIrecvScalar(comm,tag_v,nrecvs,rnodes,rlengths,&rvalues,&recv_reqs_v);
1717: }
1718: /* we are sending ni*rlengths[i] indices in each message (ni == 1 or 2, is the number of index arrays being sent) */
1719: if(ni == 2) {
1720: for (i=0; i<nrecvs; ++i) rlengths[i] *=2;
1721: }
1722: PetscCommGetNewTag(layout->comm, &tag_i);
1723: PetscPostIrecvInt(comm,tag_i,nrecvs,rnodes,rlengths,&rindices,&recv_reqs_i);
1724: PetscFree(rnodes);
1726: if(ni == 2) {
1727: for (i=0; i<nrecvs; ++i) rlengths[i] /=2;
1728: }
1730: /* prepare send buffers and offsets.
1731: sindices is the index send buffer; svalues is the value send buffer, allocated only if nv > 0.
1732: sstarts[p] gives the starting offset for values going to the pth processor, if any;
1733: because the indices (I & J, if both are being sent) are sent from the same buffer,
1734: the index offset is ni*sstarts[p].
1735: */
1736: PetscMalloc((size+1)*sizeof(PetscMPIInt),&sstarts);
1737: PetscMalloc(ni*chain->length*sizeof(PetscInt),&sindices);
1738: if(nv) {
1739: PetscMalloc(chain->length*sizeof(PetscScalar),&svalues);
1740: }
1742: /* Compute buffer offsets for the segments of data going to different processors,
1743: and zero out plengths: they will be used below as running counts when packing data
1744: into send buffers; as a result of that, plengths are recomputed by the end of the loop.
1745: */
1746: sstarts[0] = 0;
1747: for (p=0; p<size; ++p) {
1748: sstarts[p+1] = sstarts[p] + plengths[p];
1749: plengths[p] = 0;
1750: }
1752: /* Now pack the indices and, possibly, values into the appropriate buffer segments. */
1753: link = chain->first;
1754: count = 0;
1755: while(link){
1756: hunk = link->hunk;
1757: if(mask == SA_I) {
1758: ixidx = hunk->i;
1759: iyidx = hunk->j;
1760: }
1761: else {
1762: ixidx = hunk->j;
1763: iyidx = hunk->i;
1764: }
1765: val = hunk->w;
1766: for(i = 0; i < hunk->length; ++i) {
1767: p = owner[count];
1768: sindices[ni*sstarts[p]+plengths[p]] = ixidx[i];
1769: if(ni==2)
1770: sindices[ni*sstarts[p]+(sstarts[p+1]-sstarts[p])+plengths[p]] = iyidx[i];
1771: if(nv)
1772: svalues[sstarts[p]+plengths[p]] = val[i];
1773: ++plengths[p];
1774: ++count;
1775: }
1776: link = link->next;
1777: }
1778: /* Allocate send requests: for the indices, and possibly one more for the scalar values, hence +nv */
1779: PetscMalloc((1+nv)*nsends*sizeof(MPI_Request),&send_reqs);
1781: /* Post sends */
1782: for (p=0,count=0; p<size; ++p) {
1783: if (plengths[p]) {
1784: MPI_Isend(sindices+ni*sstarts[p],ni*plengths[p],MPIU_INT,p,tag_i,comm,send_reqs+count++);
1785: if(nv) {
1786: MPI_Isend(svalues+sstarts[p],plengths[p],MPIU_SCALAR,p,tag_v,comm,send_reqs+count++);
1787: }
1788: }
1789: }
1790: PetscFree2(plengths,owner);
1791: PetscFree(sstarts);
1793: /* Prepare to receive indices and values. */
1794: /* Compute the offsets of the individual received segments in the unified index/value arrays. */
1795: PetscMalloc(sizeof(PetscMPIInt)*(nrecvs+1), &rstarts);
1796: rstarts[0] = 0;
1797: for(j = 0; j < nrecvs; ++j) rstarts[j+1] = rstarts[j] + rlengths[j];
1799: alength = rstarts[nrecvs];
1800: /* Clear the SA that will store the received data segments */
1801: SAClear(achain);
1802: /* Get a hunk to pack the data into. */
1803: SAGetHunk(achain, alength, &ahunk);
1804:
1805: /* Use ahunk's data arrays as receive buffers. */
1806: if(mask == SA_I) {
1807: aixidx = ahunk->i;
1808: aiyidx = ahunk->j;
1809: }
1810: else {
1811: aiyidx = ahunk->i;
1812: aixidx = ahunk->j;
1813: }
1814: aval = ahunk->w;
1816: /* Receive indices and values, and pack them into unified arrays. */
1817: if(nv) {
1818: /* wait on scalar values receives and pack the received values into aval */
1819: count = nrecvs;
1820: rlengthtotal = 0;
1821: while (count) {
1822: PetscMPIInt n,k;
1823: MPI_Waitany(nrecvs,recv_reqs_v,&k,&recv_status);
1824: MPI_Get_count(&recv_status,MPIU_SCALAR,&n);
1825: rlengthtotal += n;
1826: count--;
1827: PetscMemcpy(aval+rstarts[k],rvalues[k],sizeof(PetscScalar)*rlengths[k]);
1828: }
1829: if (rstarts[nrecvs] != rlengthtotal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Total message lengths %D not as expected %D",rlengthtotal,rstarts[nrecvs]);
1830: }
1831:
1832: /* wait on index receives and pack the received indices into aixidx and aiyidx, as necessary. */
1833: count = nrecvs;
1834: rlengthtotal = 0;
1835: while (count) {
1836: PetscMPIInt n,k;
1837: MPI_Waitany(nrecvs,recv_reqs_i,&k,&recv_status);
1838: MPI_Get_count(&recv_status,MPIU_INT,&n);
1839: rlengthtotal += n/ni;
1840: count--;
1841: PetscMemcpy(aixidx+rstarts[k],rindices[k],sizeof(PetscInt)*rlengths[k]);
1842: if(ni == 2) {
1843: PetscMemcpy(aiyidx+rstarts[k],rindices[k]+rlengths[k],sizeof(PetscInt)*rlengths[k]);
1844: }
1845: }
1846: if (rstarts[nrecvs] != rlengthtotal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Total message lengths %D not as expected %D",rlengthtotal,rstarts[nrecvs]);
1847:
1849: PetscFree(rlengths);
1850: PetscFree(rstarts);
1851: PetscFree(rindices[0]);
1852: PetscFree(rindices);
1853: if(nv) {
1854: PetscFree(rvalues[0]);
1855: PetscFree(rvalues);
1856: }
1857: /* wait on sends */
1858: if (nsends) {
1859: PetscMalloc(sizeof(MPI_Status)*nsends,&send_statuses);
1860: MPI_Waitall(nsends,send_reqs,send_statuses);
1861: PetscFree(send_statuses);
1862: }
1864: PetscFree(sindices);
1865: if(nv) {
1866: PetscFree(svalues);
1867: }
1868: return(0);
1869: }/* SAAssemble() */