Actual source code: pack.c
1: #define PETSCDM_DLL
2:
3: #include petscda.h
4: #include petscmat.h
6: /*
7: rstart is where an array/subvector starts in the global parallel vector, so arrays
8: rstarts are meaningless (and set to the previous one) except on processor 0
9: */
11: typedef enum {VECPACK_ARRAY, VECPACK_DA, VECPACK_VECSCATTER} VecPackLinkType;
13: struct VecPackLink {
14: DA da;
15: PetscInt n,rstart; /* rstart is relative to this processor */
16: VecPackLinkType type;
17: struct VecPackLink *next;
18: };
20: typedef struct _VecPackOps *VecPackOps;
21: struct _VecPackOps {
22: PetscErrorCode (*view)(VecPack,PetscViewer);
23: PetscErrorCode (*createglobalvector)(VecPack,Vec*);
24: PetscErrorCode (*getcoloring)(VecPack,ISColoringType,ISColoring*);
25: PetscErrorCode (*getmatrix)(VecPack,MatType,Mat*);
26: PetscErrorCode (*getinterpolation)(VecPack,VecPack,Mat*,Vec*);
27: PetscErrorCode (*refine)(VecPack,MPI_Comm,VecPack*);
28: };
30: struct _p_VecPack {
31: PETSCHEADER(struct _VecPackOps);
32: PetscMPIInt rank;
33: PetscInt n,N,rstart; /* rstart is relative to all processors */
34: Vec globalvector;
35: PetscInt nDA,nredundant;
36: struct VecPackLink *next;
37: };
41: /*@C
42: VecPackCreate - Creates a vector packer, used to generate "composite"
43: vectors made up of several subvectors.
45: Collective on MPI_Comm
47: Input Parameter:
48: . comm - the processors that will share the global vector
50: Output Parameters:
51: . packer - the packer object
53: Level: advanced
55: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackScatter(),
56: VecPackGather(), VecPackCreateGlobalVector(), VecPackGetGlobalIndices(), VecPackGetAccess()
57: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
59: @*/
60: PetscErrorCode VecPackCreate(MPI_Comm comm,VecPack *packer)
61: {
63: VecPack p;
67: *packer = PETSC_NULL;
68: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
69: DMInitializePackage(PETSC_NULL);
70: #endif
72: PetscHeaderCreate(p,_p_VecPack,struct _VecPackOps,DA_COOKIE,0,"VecPack",comm,VecPackDestroy,0);
73: p->n = 0;
74: p->next = PETSC_NULL;
75: p->comm = comm;
76: p->globalvector = PETSC_NULL;
77: p->nredundant = 0;
78: p->nDA = 0;
79: MPI_Comm_rank(comm,&p->rank);
81: p->ops->createglobalvector = VecPackCreateGlobalVector;
82: p->ops->refine = VecPackRefine;
83: p->ops->getinterpolation = VecPackGetInterpolation;
84: *packer = p;
85: return(0);
86: }
90: /*@C
91: VecPackDestroy - Destroys a vector packer.
93: Collective on VecPack
95: Input Parameter:
96: . packer - the packer object
98: Level: advanced
100: .seealso VecPackCreate(), VecPackAddArray(), VecPackAddDA(), VecPackScatter(),
101: VecPackGather(), VecPackCreateGlobalVector(), VecPackGetGlobalIndices(), VecPackGetAccess()
103: @*/
104: PetscErrorCode VecPackDestroy(VecPack packer)
105: {
106: PetscErrorCode ierr;
107: struct VecPackLink *next = packer->next,*prev;
110: if (--packer->refct > 0) return(0);
111: while (next) {
112: prev = next;
113: next = next->next;
114: if (prev->type == VECPACK_DA) {
115: DADestroy(prev->da);
116: }
117: PetscFree(prev);
118: }
119: if (packer->globalvector) {
120: VecDestroy(packer->globalvector);
121: }
122: PetscHeaderDestroy(packer);
123: return(0);
124: }
126: /* --------------------------------------------------------------------------------------*/
130: PetscErrorCode VecPackGetAccess_Array(VecPack packer,struct VecPackLink *mine,Vec vec,PetscScalar **array)
131: {
133: PetscScalar *varray;
136: if (array) {
137: if (!packer->rank) {
138: VecGetArray(vec,&varray);
139: *array = varray + mine->rstart;
140: VecRestoreArray(vec,&varray);
141: } else {
142: *array = 0;
143: }
144: }
145: return(0);
146: }
150: PetscErrorCode VecPackGetAccess_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec *global)
151: {
153: PetscScalar *array;
156: if (global) {
157: DAGetGlobalVector(mine->da,global);
158: VecGetArray(vec,&array);
159: VecPlaceArray(*global,array+mine->rstart);
160: VecRestoreArray(vec,&array);
161: }
162: return(0);
163: }
167: PetscErrorCode VecPackRestoreAccess_Array(VecPack packer,struct VecPackLink *mine,Vec vec,PetscScalar **array)
168: {
170: return(0);
171: }
175: PetscErrorCode VecPackRestoreAccess_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec *global)
176: {
180: if (global) {
181: VecResetArray(*global);
182: DARestoreGlobalVector(mine->da,global);
183: }
184: return(0);
185: }
189: PetscErrorCode VecPackScatter_Array(VecPack packer,struct VecPackLink *mine,Vec vec,PetscScalar *array)
190: {
192: PetscScalar *varray;
196: if (!packer->rank) {
197: VecGetArray(vec,&varray);
198: PetscMemcpy(array,varray+mine->rstart,mine->n*sizeof(PetscScalar));
199: VecRestoreArray(vec,&varray);
200: }
201: MPI_Bcast(array,mine->n,MPIU_SCALAR,0,packer->comm);
202: return(0);
203: }
207: PetscErrorCode VecPackScatter_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec local)
208: {
210: PetscScalar *array;
211: Vec global;
214: DAGetGlobalVector(mine->da,&global);
215: VecGetArray(vec,&array);
216: VecPlaceArray(global,array+mine->rstart);
217: DAGlobalToLocalBegin(mine->da,global,INSERT_VALUES,local);
218: DAGlobalToLocalEnd(mine->da,global,INSERT_VALUES,local);
219: VecRestoreArray(vec,&array);
220: VecResetArray(global);
221: DARestoreGlobalVector(mine->da,&global);
222: return(0);
223: }
227: PetscErrorCode VecPackGather_Array(VecPack packer,struct VecPackLink *mine,Vec vec,PetscScalar *array)
228: {
230: PetscScalar *varray;
233: if (!packer->rank) {
234: VecGetArray(vec,&varray);
235: if (varray+mine->rstart == array) SETERRQ(PETSC_ERR_ARG_WRONG,"You need not VecPackGather() into objects obtained via VecPackGetAccess()");
236: PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));
237: VecRestoreArray(vec,&varray);
238: }
239: return(0);
240: }
244: PetscErrorCode VecPackGather_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec local)
245: {
247: PetscScalar *array;
248: Vec global;
251: DAGetGlobalVector(mine->da,&global);
252: VecGetArray(vec,&array);
253: VecPlaceArray(global,array+mine->rstart);
254: DALocalToGlobal(mine->da,local,INSERT_VALUES,global);
255: VecRestoreArray(vec,&array);
256: VecResetArray(global);
257: DARestoreGlobalVector(mine->da,&global);
258: return(0);
259: }
261: /* ----------------------------------------------------------------------------------*/
263: #include <stdarg.h>
267: /*@C
268: VecPackGetAccess - Allows one to access the individual packed vectors in their global
269: representation.
271: Collective on VecPack
273: Input Parameter:
274: + packer - the packer object
275: . gvec - the global vector
276: - ... - the individual sequential or parallel objects (arrays or vectors)
277:
278: Level: advanced
280: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
281: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackScatter(),
282: VecPackRestoreAccess(), VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
284: @*/
285: PetscErrorCode VecPackGetAccess(VecPack packer,Vec gvec,...)
286: {
287: va_list Argp;
288: PetscErrorCode ierr;
289: struct VecPackLink *next = packer->next;
292: if (!packer->globalvector) {
293: SETERRQ(PETSC_ERR_ORDER,"Must first create global vector with VecPackCreateGlobalVector()");
294: }
296: /* loop over packed objects, handling one at at time */
297: va_start(Argp,gvec);
298: while (next) {
299: if (next->type == VECPACK_ARRAY) {
300: PetscScalar **array;
301: array = va_arg(Argp, PetscScalar**);
302: VecPackGetAccess_Array(packer,next,gvec,array);
303: } else if (next->type == VECPACK_DA) {
304: Vec *vec;
305: vec = va_arg(Argp, Vec*);
306: VecPackGetAccess_DA(packer,next,gvec,vec);
307: } else {
308: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
309: }
310: next = next->next;
311: }
312: va_end(Argp);
313: return(0);
314: }
318: /*@C
319: VecPackRestoreAccess - Allows one to access the individual packed vectors in their global
320: representation.
322: Collective on VecPack
324: Input Parameter:
325: + packer - the packer object
326: . gvec - the global vector
327: - ... - the individual sequential or parallel objects (arrays or vectors)
328:
329: Level: advanced
331: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
332: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackScatter(),
333: VecPackRestoreAccess()
335: @*/
336: PetscErrorCode VecPackRestoreAccess(VecPack packer,Vec gvec,...)
337: {
338: va_list Argp;
339: PetscErrorCode ierr;
340: struct VecPackLink *next = packer->next;
343: if (!packer->globalvector) {
344: SETERRQ(PETSC_ERR_ORDER,"Must first create global vector with VecPackCreateGlobalVector()");
345: }
347: /* loop over packed objects, handling one at at time */
348: va_start(Argp,gvec);
349: while (next) {
350: if (next->type == VECPACK_ARRAY) {
351: PetscScalar **array;
352: array = va_arg(Argp, PetscScalar**);
353: VecPackRestoreAccess_Array(packer,next,gvec,array);
354: } else if (next->type == VECPACK_DA) {
355: Vec *vec;
356: vec = va_arg(Argp, Vec*);
357: VecPackRestoreAccess_DA(packer,next,gvec,vec);
358: } else {
359: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
360: }
361: next = next->next;
362: }
363: va_end(Argp);
364: return(0);
365: }
369: /*@C
370: VecPackScatter - Scatters from a global packed vector into its individual local vectors
372: Collective on VecPack
374: Input Parameter:
375: + packer - the packer object
376: . gvec - the global vector
377: - ... - the individual sequential objects (arrays or vectors)
378:
379: Level: advanced
381: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
382: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
383: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
385: @*/
386: PetscErrorCode VecPackScatter(VecPack packer,Vec gvec,...)
387: {
388: va_list Argp;
389: PetscErrorCode ierr;
390: struct VecPackLink *next = packer->next;
393: if (!packer->globalvector) {
394: SETERRQ(PETSC_ERR_ORDER,"Must first create global vector with VecPackCreateGlobalVector()");
395: }
397: /* loop over packed objects, handling one at at time */
398: va_start(Argp,gvec);
399: while (next) {
400: if (next->type == VECPACK_ARRAY) {
401: PetscScalar *array;
402: array = va_arg(Argp, PetscScalar*);
403: VecPackScatter_Array(packer,next,gvec,array);
404: } else if (next->type == VECPACK_DA) {
405: Vec vec;
406: vec = va_arg(Argp, Vec);
408: VecPackScatter_DA(packer,next,gvec,vec);
409: } else {
410: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
411: }
412: next = next->next;
413: }
414: va_end(Argp);
415: return(0);
416: }
420: /*@C
421: VecPackGather - Gathers into a global packed vector from its individual local vectors
423: Collective on VecPack
425: Input Parameter:
426: + packer - the packer object
427: . gvec - the global vector
428: - ... - the individual sequential objects (arrays or vectors)
429:
430: Level: advanced
432: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
433: VecPackScatter(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
434: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
436: @*/
437: PetscErrorCode VecPackGather(VecPack packer,Vec gvec,...)
438: {
439: va_list Argp;
440: PetscErrorCode ierr;
441: struct VecPackLink *next = packer->next;
444: if (!packer->globalvector) {
445: SETERRQ(PETSC_ERR_ORDER,"Must first create global vector with VecPackCreateGlobalVector()");
446: }
448: /* loop over packed objects, handling one at at time */
449: va_start(Argp,gvec);
450: while (next) {
451: if (next->type == VECPACK_ARRAY) {
452: PetscScalar *array;
453: array = va_arg(Argp, PetscScalar*);
454: VecPackGather_Array(packer,next,gvec,array);
455: } else if (next->type == VECPACK_DA) {
456: Vec vec;
457: vec = va_arg(Argp, Vec);
459: VecPackGather_DA(packer,next,gvec,vec);
460: } else {
461: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
462: }
463: next = next->next;
464: }
465: va_end(Argp);
466: return(0);
467: }
471: /*@C
472: VecPackAddArray - adds an "redundant" array to a VecPack. The array values will
473: be stored in part of the array on processor 0.
475: Collective on VecPack
477: Input Parameter:
478: + packer - the packer object
479: - n - the length of the array
480:
481: Level: advanced
483: .seealso VecPackDestroy(), VecPackGather(), VecPackAddDA(), VecPackCreateGlobalVector(),
484: VecPackScatter(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
485: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
487: @*/
488: PetscErrorCode VecPackAddArray(VecPack packer,PetscInt n)
489: {
490: struct VecPackLink *mine,*next = packer->next;
491: PetscErrorCode ierr;
494: if (packer->globalvector) {
495: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot add an array once you have called VecPackCreateGlobalVector()");
496: }
498: /* create new link */
499: PetscNew(struct VecPackLink,&mine);
500: mine->n = n;
501: mine->da = PETSC_NULL;
502: mine->type = VECPACK_ARRAY;
503: mine->next = PETSC_NULL;
504: if (!packer->rank) packer->n += n;
506: /* add to end of list */
507: if (!next) {
508: packer->next = mine;
509: } else {
510: while (next->next) next = next->next;
511: next->next = mine;
512: }
513: packer->nredundant++;
514: return(0);
515: }
519: /*@C
520: VecPackAddDA - adds a DA vector to a VecPack
522: Collective on VecPack
524: Input Parameter:
525: + packer - the packer object
526: - da - the DA object
527:
528: Level: advanced
530: .seealso VecPackDestroy(), VecPackGather(), VecPackAddDA(), VecPackCreateGlobalVector(),
531: VecPackScatter(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
532: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
534: @*/
535: PetscErrorCode VecPackAddDA(VecPack packer,DA da)
536: {
537: PetscErrorCode ierr;
538: PetscInt n;
539: struct VecPackLink *mine,*next = packer->next;
540: Vec global;
543: if (packer->globalvector) {
544: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DA once you have called VecPackCreateGlobalVector()");
545: }
547: /* create new link */
548: PetscNew(struct VecPackLink,&mine);
549: PetscObjectReference((PetscObject)da);
550: DAGetGlobalVector(da,&global);
551: VecGetLocalSize(global,&n);
552: DARestoreGlobalVector(da,&global);
553: mine->n = n;
554: mine->da = da;
555: mine->type = VECPACK_DA;
556: mine->next = PETSC_NULL;
557: packer->n += n;
559: /* add to end of list */
560: if (!next) {
561: packer->next = mine;
562: } else {
563: while (next->next) next = next->next;
564: next->next = mine;
565: }
566: packer->nDA++;
567: return(0);
568: }
572: /*@C
573: VecPackCreateGlobalVector - Creates a vector of the correct size to be gathered into
574: by the packer.
576: Collective on VecPack
578: Input Parameter:
579: . packer - the packer object
581: Output Parameters:
582: . gvec - the global vector
584: Level: advanced
586: Notes: Once this has been created you cannot add additional arrays or vectors to be packed.
588: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackScatter(),
589: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
590: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
592: @*/
593: PetscErrorCode VecPackCreateGlobalVector(VecPack packer,Vec *gvec)
594: {
595: PetscErrorCode ierr;
596: PetscInt nprev = 0;
597: PetscMPIInt rank;
598: struct VecPackLink *next = packer->next;
601: if (packer->globalvector) {
602: VecDuplicate(packer->globalvector,gvec);
603: } else {
604: VecCreateMPI(packer->comm,packer->n,PETSC_DETERMINE,gvec);
605: PetscObjectReference((PetscObject)*gvec);
606: packer->globalvector = *gvec;
608: VecGetSize(*gvec,&packer->N);
609: VecGetOwnershipRange(*gvec,&packer->rstart,PETSC_NULL);
610:
611: /* now set the rstart for each linked array/vector */
612: MPI_Comm_rank(packer->comm,&rank);
613: while (next) {
614: next->rstart = nprev;
615: if (!rank || next->type != VECPACK_ARRAY) nprev += next->n;
616: next = next->next;
617: }
618: }
619: return(0);
620: }
624: /*@C
625: VecPackGetGlobalIndices - Gets the global indices for all the entries in the packed
626: vectors.
628: Collective on VecPack
630: Input Parameter:
631: . packer - the packer object
633: Output Parameters:
634: . idx - the individual indices for each packed vector/array. Note that this includes
635: all the ghost points that individual ghosted DA's may have.
636:
637: Level: advanced
639: Notes:
640: The idx parameters should be freed by the calling routine with PetscFree()
642: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
643: VecPackGather(), VecPackCreate(), VecPackGetAccess(),
644: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
646: @*/
647: PetscErrorCode VecPackGetGlobalIndices(VecPack packer,...)
648: {
649: va_list Argp;
650: PetscErrorCode ierr;
651: PetscInt i,**idx,n;
652: struct VecPackLink *next = packer->next;
653: Vec global,dglobal;
654: PF pf;
655: PetscScalar *array;
658: VecPackCreateGlobalVector(packer,&global);
660: /* put 0 to N-1 into the global vector */
661: PFCreate(PETSC_COMM_WORLD,1,1,&pf);
662: PFSetType(pf,PFIDENTITY,PETSC_NULL);
663: PFApplyVec(pf,PETSC_NULL,global);
664: PFDestroy(pf);
666: /* loop over packed objects, handling one at at time */
667: va_start(Argp,packer);
668: while (next) {
669: idx = va_arg(Argp, PetscInt**);
671: if (next->type == VECPACK_ARRAY) {
672:
673: PetscMalloc(next->n*sizeof(PetscInt),idx);
674: if (!packer->rank) {
675: VecGetArray(global,&array);
676: array += next->rstart;
677: for (i=0; i<next->n; i++) (*idx)[i] = (PetscInt)PetscRealPart(array[i]);
678: array -= next->rstart;
679: VecRestoreArray(global,&array);
680: }
681: MPI_Bcast(*idx,next->n,MPIU_INT,0,packer->comm);
683: } else if (next->type == VECPACK_DA) {
684: Vec local;
686: DACreateLocalVector(next->da,&local);
687: VecGetArray(global,&array);
688: array += next->rstart;
689: DAGetGlobalVector(next->da,&dglobal);
690: VecPlaceArray(dglobal,array);
691: DAGlobalToLocalBegin(next->da,dglobal,INSERT_VALUES,local);
692: DAGlobalToLocalEnd(next->da,dglobal,INSERT_VALUES,local);
693: array -= next->rstart;
694: VecRestoreArray(global,&array);
695: VecResetArray(dglobal);
696: DARestoreGlobalVector(next->da,&dglobal);
698: VecGetArray(local,&array);
699: VecGetSize(local,&n);
700: PetscMalloc(n*sizeof(PetscInt),idx);
701: for (i=0; i<n; i++) (*idx)[i] = (PetscInt)PetscRealPart(array[i]);
702: VecRestoreArray(local,&array);
703: VecDestroy(local);
705: } else {
706: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
707: }
708: next = next->next;
709: }
710: va_end(Argp);
711: VecDestroy(global);
712: return(0);
713: }
715: /* -------------------------------------------------------------------------------------*/
718: PetscErrorCode VecPackGetLocalVectors_Array(VecPack packer,struct VecPackLink *mine,PetscScalar **array)
719: {
722: PetscMalloc(mine->n*sizeof(PetscScalar),array);
723: return(0);
724: }
728: PetscErrorCode VecPackGetLocalVectors_DA(VecPack packer,struct VecPackLink *mine,Vec *local)
729: {
732: DAGetLocalVector(mine->da,local);
733: return(0);
734: }
738: PetscErrorCode VecPackRestoreLocalVectors_Array(VecPack packer,struct VecPackLink *mine,PetscScalar **array)
739: {
742: PetscFree(*array);
743: return(0);
744: }
748: PetscErrorCode VecPackRestoreLocalVectors_DA(VecPack packer,struct VecPackLink *mine,Vec *local)
749: {
752: DARestoreLocalVector(mine->da,local);
753: return(0);
754: }
758: /*@C
759: VecPackGetLocalVectors - Gets local vectors and arrays for each part of a VecPack.'
760: Use VecPakcRestoreLocalVectors() to return them.
762: Collective on VecPack
764: Input Parameter:
765: . packer - the packer object
766:
767: Output Parameter:
768: . ... - the individual sequential objects (arrays or vectors)
769:
770: Level: advanced
772: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
773: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
774: VecPackRestoreLocalVectors()
776: @*/
777: PetscErrorCode VecPackGetLocalVectors(VecPack packer,...)
778: {
779: va_list Argp;
780: PetscErrorCode ierr;
781: struct VecPackLink *next = packer->next;
785: /* loop over packed objects, handling one at at time */
786: va_start(Argp,packer);
787: while (next) {
788: if (next->type == VECPACK_ARRAY) {
789: PetscScalar **array;
790: array = va_arg(Argp, PetscScalar**);
791: VecPackGetLocalVectors_Array(packer,next,array);
792: } else if (next->type == VECPACK_DA) {
793: Vec *vec;
794: vec = va_arg(Argp, Vec*);
795: VecPackGetLocalVectors_DA(packer,next,vec);
796: } else {
797: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
798: }
799: next = next->next;
800: }
801: va_end(Argp);
802: return(0);
803: }
807: /*@C
808: VecPackRestoreLocalVectors - Restores local vectors and arrays for each part of a VecPack.'
809: Use VecPakcRestoreLocalVectors() to return them.
811: Collective on VecPack
813: Input Parameter:
814: . packer - the packer object
815:
816: Output Parameter:
817: . ... - the individual sequential objects (arrays or vectors)
818:
819: Level: advanced
821: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
822: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
823: VecPackGetLocalVectors()
825: @*/
826: PetscErrorCode VecPackRestoreLocalVectors(VecPack packer,...)
827: {
828: va_list Argp;
829: PetscErrorCode ierr;
830: struct VecPackLink *next = packer->next;
834: /* loop over packed objects, handling one at at time */
835: va_start(Argp,packer);
836: while (next) {
837: if (next->type == VECPACK_ARRAY) {
838: PetscScalar **array;
839: array = va_arg(Argp, PetscScalar**);
840: VecPackRestoreLocalVectors_Array(packer,next,array);
841: } else if (next->type == VECPACK_DA) {
842: Vec *vec;
843: vec = va_arg(Argp, Vec*);
844: VecPackRestoreLocalVectors_DA(packer,next,vec);
845: } else {
846: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
847: }
848: next = next->next;
849: }
850: va_end(Argp);
851: return(0);
852: }
854: /* -------------------------------------------------------------------------------------*/
857: PetscErrorCode VecPackGetEntries_Array(VecPack packer,struct VecPackLink *mine,PetscInt *n)
858: {
860: if (n) *n = mine->n;
861: return(0);
862: }
866: PetscErrorCode VecPackGetEntries_DA(VecPack packer,struct VecPackLink *mine,DA *da)
867: {
869: if (da) *da = mine->da;
870: return(0);
871: }
875: /*@C
876: VecPackGetEntries - Gets the DA, redundant size, etc for each entry in a VecPack.
877: Use VecPackRestoreEntries() to return them.
879: Collective on VecPack
881: Input Parameter:
882: . packer - the packer object
883:
884: Output Parameter:
885: . ... - the individual entries, DAs or integer sizes)
886:
887: Level: advanced
889: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
890: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
891: VecPackRestoreLocalVectors(), VecPackGetLocalVectors(), VecPackRestoreEntries(),
892: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
894: @*/
895: PetscErrorCode VecPackGetEntries(VecPack packer,...)
896: {
897: va_list Argp;
898: PetscErrorCode ierr;
899: struct VecPackLink *next = packer->next;
903: /* loop over packed objects, handling one at at time */
904: va_start(Argp,packer);
905: while (next) {
906: if (next->type == VECPACK_ARRAY) {
907: PetscInt *n;
908: n = va_arg(Argp, PetscInt*);
909: VecPackGetEntries_Array(packer,next,n);
910: } else if (next->type == VECPACK_DA) {
911: DA *da;
912: da = va_arg(Argp, DA*);
913: VecPackGetEntries_DA(packer,next,da);
914: } else {
915: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
916: }
917: next = next->next;
918: }
919: va_end(Argp);
920: return(0);
921: }
925: /*@C
926: VecPackRefine - Refines a VecPack by refining all of its DAs
928: Collective on VecPack
930: Input Parameters:
931: + packer - the packer object
932: - comm - communicator to contain the new DM object, usually PETSC_NULL
934: Output Parameter:
935: . fine - new packer
936:
937: Level: advanced
939: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
940: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
941: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
943: @*/
944: PetscErrorCode VecPackRefine(VecPack packer,MPI_Comm comm,VecPack *fine)
945: {
946: PetscErrorCode ierr;
947: struct VecPackLink *next = packer->next;
948: DA da;
951: VecPackCreate(comm,fine);
953: /* loop over packed objects, handling one at at time */
954: while (next) {
955: if (next->type == VECPACK_ARRAY) {
956: VecPackAddArray(*fine,next->n);
957: } else if (next->type == VECPACK_DA) {
958: DARefine(next->da,comm,&da);
959: VecPackAddDA(*fine,da);
960: PetscObjectDereference((PetscObject)da);
961: } else {
962: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
963: }
964: next = next->next;
965: }
966: return(0);
967: }
969: #include petscmat.h
971: struct MatPackLink {
972: Mat A;
973: struct MatPackLink *next;
974: };
976: struct MatPack {
977: VecPack right,left;
978: struct MatPackLink *next;
979: };
983: PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscTruth add)
984: {
985: struct MatPack *mpack;
986: struct VecPackLink *xnext,*ynext;
987: struct MatPackLink *anext;
988: PetscScalar *xarray,*yarray;
989: PetscErrorCode ierr;
990: PetscInt i;
991: Vec xglobal,yglobal;
994: MatShellGetContext(A,(void**)&mpack);
995: xnext = mpack->right->next;
996: ynext = mpack->left->next;
997: anext = mpack->next;
999: while (xnext) {
1000: if (xnext->type == VECPACK_ARRAY) {
1001: if (!mpack->right->rank) {
1002: VecGetArray(x,&xarray);
1003: VecGetArray(y,&yarray);
1004: if (add) {
1005: for (i=0; i<xnext->n; i++) {
1006: yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
1007: }
1008: } else {
1009: PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));
1010: }
1011: VecRestoreArray(x,&xarray);
1012: VecRestoreArray(y,&yarray);
1013: }
1014: } else if (xnext->type == VECPACK_DA) {
1015: VecGetArray(x,&xarray);
1016: VecGetArray(y,&yarray);
1017: DAGetGlobalVector(xnext->da,&xglobal);
1018: DAGetGlobalVector(ynext->da,&yglobal);
1019: VecPlaceArray(xglobal,xarray+xnext->rstart);
1020: VecPlaceArray(yglobal,yarray+ynext->rstart);
1021: if (add) {
1022: MatMultAdd(anext->A,xglobal,yglobal,yglobal);
1023: } else {
1024: MatMult(anext->A,xglobal,yglobal);
1025: }
1026: VecRestoreArray(x,&xarray);
1027: VecRestoreArray(y,&yarray);
1028: VecResetArray(xglobal);
1029: VecResetArray(yglobal);
1030: DARestoreGlobalVector(xnext->da,&xglobal);
1031: DARestoreGlobalVector(ynext->da,&yglobal);
1032: anext = anext->next;
1033: } else {
1034: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
1035: }
1036: xnext = xnext->next;
1037: ynext = ynext->next;
1038: }
1039: return(0);
1040: }
1044: PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
1045: {
1048: if (z != y) SETERRQ(PETSC_ERR_SUP,"Handles y == z only");
1049: MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);
1050: return(0);
1051: }
1055: PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y)
1056: {
1059: MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);
1060: return(0);
1061: }
1065: PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
1066: {
1067: struct MatPack *mpack;
1068: struct VecPackLink *xnext,*ynext;
1069: struct MatPackLink *anext;
1070: PetscScalar *xarray,*yarray;
1071: PetscErrorCode ierr;
1072: Vec xglobal,yglobal;
1075: MatShellGetContext(A,(void**)&mpack);
1076: xnext = mpack->left->next;
1077: ynext = mpack->right->next;
1078: anext = mpack->next;
1080: while (xnext) {
1081: if (xnext->type == VECPACK_ARRAY) {
1082: if (!mpack->right->rank) {
1083: VecGetArray(x,&xarray);
1084: VecGetArray(y,&yarray);
1085: PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));
1086: VecRestoreArray(x,&xarray);
1087: VecRestoreArray(y,&yarray);
1088: }
1089: } else if (xnext->type == VECPACK_DA) {
1090: VecGetArray(x,&xarray);
1091: VecGetArray(y,&yarray);
1092: DAGetGlobalVector(xnext->da,&xglobal);
1093: DAGetGlobalVector(ynext->da,&yglobal);
1094: VecPlaceArray(xglobal,xarray+xnext->rstart);
1095: VecPlaceArray(yglobal,yarray+ynext->rstart);
1096: MatMultTranspose(anext->A,xglobal,yglobal);
1097: VecRestoreArray(x,&xarray);
1098: VecRestoreArray(y,&yarray);
1099: VecResetArray(xglobal);
1100: VecResetArray(yglobal);
1101: DARestoreGlobalVector(xnext->da,&xglobal);
1102: DARestoreGlobalVector(ynext->da,&yglobal);
1103: anext = anext->next;
1104: } else {
1105: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
1106: }
1107: xnext = xnext->next;
1108: ynext = ynext->next;
1109: }
1110: return(0);
1111: }
1115: PetscErrorCode MatDestroy_Shell_Pack(Mat A)
1116: {
1117: struct MatPack *mpack;
1118: struct MatPackLink *anext,*oldanext;
1119: PetscErrorCode ierr;
1122: MatShellGetContext(A,(void**)&mpack);
1123: anext = mpack->next;
1125: while (anext) {
1126: MatDestroy(anext->A);
1127: oldanext = anext;
1128: anext = anext->next;
1129: PetscFree(oldanext);
1130: }
1131: PetscFree(mpack);
1132: PetscObjectChangeTypeName((PetscObject)A,0);
1133: return(0);
1134: }
1138: /*@C
1139: VecPackGetInterpolation - GetInterpolations a VecPack by refining all of its DAs
1141: Collective on VecPack
1143: Input Parameters:
1144: + coarse - coarse grid packer
1145: - fine - fine grid packer
1147: Output Parameter:
1148: + A - interpolation matrix
1149: - v - scaling vector
1150:
1151: Level: advanced
1153: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
1154: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
1155: VecPackGetLocalVectors(), VecPackRestoreLocalVectors()
1157: @*/
1158: PetscErrorCode VecPackGetInterpolation(VecPack coarse,VecPack fine,Mat *A,Vec *v)
1159: {
1160: PetscErrorCode ierr;
1161: PetscInt m,n,M,N;
1162: struct VecPackLink *nextc = coarse->next;
1163: struct VecPackLink *nextf = fine->next;
1164: struct MatPackLink *nextmat,*pnextmat = 0;
1165: struct MatPack *mpack;
1166: Vec gcoarse,gfine;
1169: /* use global vectors only for determining matrix layout */
1170: VecPackCreateGlobalVector(coarse,&gcoarse);
1171: VecPackCreateGlobalVector(fine,&gfine);
1172: VecGetLocalSize(gcoarse,&n);
1173: VecGetLocalSize(gfine,&m);
1174: VecGetSize(gcoarse,&N);
1175: VecGetSize(gfine,&M);
1176: VecDestroy(gcoarse);
1177: VecDestroy(gfine);
1179: PetscNew(struct MatPack,&mpack);
1180: mpack->right = coarse;
1181: mpack->left = fine;
1182: MatCreate(fine->comm,A);
1183: MatSetSizes(*A,m,n,M,N);
1184: MatSetType(*A,MATSHELL);
1185: MatShellSetContext(*A,mpack);
1186: MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);
1187: MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);
1188: MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);
1189: MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);
1191: /* loop over packed objects, handling one at at time */
1192: while (nextc) {
1193: if (nextc->type != nextf->type) SETERRQ(PETSC_ERR_ARG_INCOMP,"Two VecPack have different layout");
1195: if (nextc->type == VECPACK_ARRAY) {
1196: ;
1197: } else if (nextc->type == VECPACK_DA) {
1198: PetscNew(struct MatPackLink,&nextmat);
1199: nextmat->next = 0;
1200: if (pnextmat) {
1201: pnextmat->next = nextmat;
1202: pnextmat = nextmat;
1203: } else {
1204: pnextmat = nextmat;
1205: mpack->next = nextmat;
1206: }
1207: DAGetInterpolation(nextc->da,nextf->da,&nextmat->A,PETSC_NULL);
1208: } else {
1209: SETERRQ(PETSC_ERR_SUP,"Cannot handle that object type yet");
1210: }
1211: nextc = nextc->next;
1212: nextf = nextf->next;
1213: }
1214: return(0);
1215: }