Actual source code: pack.c
1: /*$Id: pack.c,v 1.19 2001/08/07 03:04:45 balay Exp $*/
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: int 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: int (*view)(VecPack,PetscViewer);
23: int (*createglobalvector)(VecPack,Vec*);
24: int (*getcoloring)(VecPack,ISColoringType,ISColoring*);
25: int (*getmatrix)(VecPack,MatType,Mat*);
26: int (*getinterpolation)(VecPack,VecPack,Mat*,Vec*);
27: int (*refine)(VecPack,MPI_Comm,VecPack*);
28: };
30: struct _p_VecPack {
31: PETSCHEADER(struct _VecPackOps)
32: int rank;
33: int n,N,rstart; /* rstart is relative to all processors */
34: Vec globalvector;
35: int nDA,nredundant;
36: struct VecPackLink *next;
37: };
39: /*@C
40: VecPackCreate - Creates a vector packer, used to generate "composite"
41: vectors made up of several subvectors.
43: Collective on MPI_Comm
45: Input Parameter:
46: . comm - the processors that will share the global vector
48: Output Parameters:
49: . packer - the packer object
51: Level: advanced
53: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackScatter(),
54: VecPackGather(), VecPackCreateGlobalVector(), VecPackGetGlobalIndices(), VecPackGetAccess()
56: @*/
57: int VecPackCreate(MPI_Comm comm,VecPack *packer)
58: {
59: int ierr;
60: VecPack p;
64: *packer = PETSC_NULL;
65: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
66: DMInitializePackage(PETSC_NULL);
67: #endif
69: PetscHeaderCreate(p,_p_VecPack,struct _VecPackOps,DA_COOKIE,0,"VecPack",comm,VecPackDestroy,0);
70: PetscLogObjectCreate(p);
71: p->n = 0;
72: p->next = PETSC_NULL;
73: p->comm = comm;
74: p->globalvector = PETSC_NULL;
75: p->nredundant = 0;
76: p->nDA = 0;
77: ierr = MPI_Comm_rank(comm,&p->rank);
79: p->ops->createglobalvector = VecPackCreateGlobalVector;
80: p->ops->refine = VecPackRefine;
81: p->ops->getinterpolation = VecPackGetInterpolation;
82: *packer = p;
83: return(0);
84: }
86: /*@C
87: VecPackDestroy - Destroys a vector packer.
89: Collective on VecPack
91: Input Parameter:
92: . packer - the packer object
94: Level: advanced
96: .seealso VecPackCreate(), VecPackAddArray(), VecPackAddDA(), VecPackScatter(),
97: VecPackGather(), VecPackCreateGlobalVector(), VecPackGetGlobalIndices(), VecPackGetAccess()
99: @*/
100: int VecPackDestroy(VecPack packer)
101: {
102: int ierr;
103: struct VecPackLink *next = packer->next,*prev;
106: if (--packer->refct > 0) return(0);
107: while (next) {
108: prev = next;
109: next = next->next;
110: if (prev->type == VECPACK_DA) {
111: DADestroy(prev->da);
112: }
113: PetscFree(prev);
114: }
115: if (packer->globalvector) {
116: VecDestroy(packer->globalvector);
117: }
118: PetscHeaderDestroy(packer);
119: return(0);
120: }
122: /* --------------------------------------------------------------------------------------*/
124: int VecPackGetAccess_Array(VecPack packer,struct VecPackLink *mine,Vec vec,PetscScalar **array)
125: {
126: int ierr;
127: PetscScalar *varray;
130: if (array) {
131: if (!packer->rank) {
132: ierr = VecGetArray(vec,&varray);
133: *array = varray + mine->rstart;
134: ierr = VecRestoreArray(vec,&varray);
135: } else {
136: *array = 0;
137: }
138: }
139: return(0);
140: }
142: int VecPackGetAccess_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec *global)
143: {
144: int ierr;
145: PetscScalar *array;
148: if (global) {
149: ierr = DAGetGlobalVector(mine->da,global);
150: ierr = VecGetArray(vec,&array);
151: ierr = VecPlaceArray(*global,array+mine->rstart);
152: ierr = VecRestoreArray(vec,&array);
153: }
154: return(0);
155: }
157: int VecPackRestoreAccess_Array(VecPack packer,struct VecPackLink *mine,Vec vec,PetscScalar **array)
158: {
160: return(0);
161: }
163: int VecPackRestoreAccess_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec *global)
164: {
165: int ierr;
168: if (global) {
169: VecResetArray(*global);
170: DARestoreGlobalVector(mine->da,global);
171: }
172: return(0);
173: }
175: int VecPackScatter_Array(VecPack packer,struct VecPackLink *mine,Vec vec,PetscScalar *array)
176: {
177: int ierr;
178: PetscScalar *varray;
182: if (!packer->rank) {
183: ierr = VecGetArray(vec,&varray);
184: ierr = PetscMemcpy(array,varray+mine->rstart,mine->n*sizeof(PetscScalar));
185: ierr = VecRestoreArray(vec,&varray);
186: }
187: ierr = MPI_Bcast(array,mine->n,MPIU_SCALAR,0,packer->comm);
188: return(0);
189: }
191: int VecPackScatter_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec local)
192: {
193: int ierr;
194: PetscScalar *array;
195: Vec global;
198: DAGetGlobalVector(mine->da,&global);
199: VecGetArray(vec,&array);
200: VecPlaceArray(global,array+mine->rstart);
201: DAGlobalToLocalBegin(mine->da,global,INSERT_VALUES,local);
202: DAGlobalToLocalEnd(mine->da,global,INSERT_VALUES,local);
203: VecRestoreArray(vec,&array);
204: VecResetArray(global);
205: DARestoreGlobalVector(mine->da,&global);
206: return(0);
207: }
209: int VecPackGather_Array(VecPack packer,struct VecPackLink *mine,Vec vec,PetscScalar *array)
210: {
211: int ierr;
212: PetscScalar *varray;
215: if (!packer->rank) {
216: ierr = VecGetArray(vec,&varray);
217: if (varray+mine->rstart == array) SETERRQ(1,"You need not VecPackGather() into objects obtained via VecPackGetAccess()");
218: ierr = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));
219: ierr = VecRestoreArray(vec,&varray);
220: }
221: return(0);
222: }
224: int VecPackGather_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec local)
225: {
226: int ierr;
227: PetscScalar *array;
228: Vec global;
231: DAGetGlobalVector(mine->da,&global);
232: VecGetArray(vec,&array);
233: VecPlaceArray(global,array+mine->rstart);
234: DALocalToGlobal(mine->da,local,INSERT_VALUES,global);
235: VecRestoreArray(vec,&array);
236: VecResetArray(global);
237: DARestoreGlobalVector(mine->da,&global);
238: return(0);
239: }
241: /* ----------------------------------------------------------------------------------*/
243: #include <stdarg.h>
245: /*@C
246: VecPackGetAccess - Allows one to access the individual packed vectors in their global
247: representation.
249: Collective on VecPack
251: Input Parameter:
252: + packer - the packer object
253: . gvec - the global vector
254: - ... - the individual sequential or parallel objects (arrays or vectors)
255:
256: Level: advanced
258: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
259: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackScatter(),
260: VecPackRestoreAccess()
262: @*/
263: int VecPackGetAccess(VecPack packer,Vec gvec,...)
264: {
265: va_list Argp;
266: int ierr;
267: struct VecPackLink *next = packer->next;
270: if (!packer->globalvector) {
271: SETERRQ(1,"Must first create global vector with VecPackCreateGlobalVector()");
272: }
274: /* loop over packed objects, handling one at at time */
275: va_start(Argp,gvec);
276: while (next) {
277: if (next->type == VECPACK_ARRAY) {
278: PetscScalar **array;
279: array = va_arg(Argp, PetscScalar**);
280: ierr = VecPackGetAccess_Array(packer,next,gvec,array);
281: } else if (next->type == VECPACK_DA) {
282: Vec *vec;
283: vec = va_arg(Argp, Vec*);
284: VecPackGetAccess_DA(packer,next,gvec,vec);
285: } else {
286: SETERRQ(1,"Cannot handle that object type yet");
287: }
288: next = next->next;
289: }
290: va_end(Argp);
291: return(0);
292: }
294: /*@C
295: VecPackRestoreAccess - Allows one to access the individual packed vectors in their global
296: representation.
298: Collective on VecPack
300: Input Parameter:
301: + packer - the packer object
302: . gvec - the global vector
303: - ... - the individual sequential or parallel objects (arrays or vectors)
304:
305: Level: advanced
307: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
308: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackScatter(),
309: VecPackRestoreAccess()
311: @*/
312: int VecPackRestoreAccess(VecPack packer,Vec gvec,...)
313: {
314: va_list Argp;
315: int ierr;
316: struct VecPackLink *next = packer->next;
319: if (!packer->globalvector) {
320: SETERRQ(1,"Must first create global vector with VecPackCreateGlobalVector()");
321: }
323: /* loop over packed objects, handling one at at time */
324: va_start(Argp,gvec);
325: while (next) {
326: if (next->type == VECPACK_ARRAY) {
327: PetscScalar **array;
328: array = va_arg(Argp, PetscScalar**);
329: ierr = VecPackRestoreAccess_Array(packer,next,gvec,array);
330: } else if (next->type == VECPACK_DA) {
331: Vec *vec;
332: vec = va_arg(Argp, Vec*);
333: VecPackRestoreAccess_DA(packer,next,gvec,vec);
334: } else {
335: SETERRQ(1,"Cannot handle that object type yet");
336: }
337: next = next->next;
338: }
339: va_end(Argp);
340: return(0);
341: }
343: /*@C
344: VecPackScatter - Scatters from a global packed vector into its individual local vectors
346: Collective on VecPack
348: Input Parameter:
349: + packer - the packer object
350: . gvec - the global vector
351: - ... - the individual sequential objects (arrays or vectors)
352:
353: Level: advanced
355: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
356: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()
358: @*/
359: int VecPackScatter(VecPack packer,Vec gvec,...)
360: {
361: va_list Argp;
362: int ierr;
363: struct VecPackLink *next = packer->next;
366: if (!packer->globalvector) {
367: SETERRQ(1,"Must first create global vector with VecPackCreateGlobalVector()");
368: }
370: /* loop over packed objects, handling one at at time */
371: va_start(Argp,gvec);
372: while (next) {
373: if (next->type == VECPACK_ARRAY) {
374: PetscScalar *array;
375: array = va_arg(Argp, PetscScalar*);
376: VecPackScatter_Array(packer,next,gvec,array);
377: } else if (next->type == VECPACK_DA) {
378: Vec vec;
379: vec = va_arg(Argp, Vec);
381: VecPackScatter_DA(packer,next,gvec,vec);
382: } else {
383: SETERRQ(1,"Cannot handle that object type yet");
384: }
385: next = next->next;
386: }
387: va_end(Argp);
388: return(0);
389: }
391: /*@C
392: VecPackGather - Gathers into a global packed vector from its individual local vectors
394: Collective on VecPack
396: Input Parameter:
397: + packer - the packer object
398: . gvec - the global vector
399: - ... - the individual sequential objects (arrays or vectors)
400:
401: Level: advanced
403: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
404: VecPackScatter(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()
406: @*/
407: int VecPackGather(VecPack packer,Vec gvec,...)
408: {
409: va_list Argp;
410: int ierr;
411: struct VecPackLink *next = packer->next;
414: if (!packer->globalvector) {
415: SETERRQ(1,"Must first create global vector with VecPackCreateGlobalVector()");
416: }
418: /* loop over packed objects, handling one at at time */
419: va_start(Argp,gvec);
420: while (next) {
421: if (next->type == VECPACK_ARRAY) {
422: PetscScalar *array;
423: array = va_arg(Argp, PetscScalar*);
424: ierr = VecPackGather_Array(packer,next,gvec,array);
425: } else if (next->type == VECPACK_DA) {
426: Vec vec;
427: vec = va_arg(Argp, Vec);
429: VecPackGather_DA(packer,next,gvec,vec);
430: } else {
431: SETERRQ(1,"Cannot handle that object type yet");
432: }
433: next = next->next;
434: }
435: va_end(Argp);
436: return(0);
437: }
439: /*@C
440: VecPackAddArray - adds an "redundant" array to a VecPack. The array values will
441: be stored in part of the array on processor 0.
443: Collective on VecPack
445: Input Parameter:
446: + packer - the packer object
447: - n - the length of the array
448:
449: Level: advanced
451: .seealso VecPackDestroy(), VecPackGather(), VecPackAddDA(), VecPackCreateGlobalVector(),
452: VecPackScatter(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()
454: @*/
455: int VecPackAddArray(VecPack packer,int n)
456: {
457: struct VecPackLink *mine,*next = packer->next;
461: if (packer->globalvector) {
462: SETERRQ(1,"Cannot add an array once you have called VecPackCreateGlobalVector()");
463: }
465: /* create new link */
466: ierr = PetscNew(struct VecPackLink,&mine);
467: mine->n = n;
468: mine->da = PETSC_NULL;
469: mine->type = VECPACK_ARRAY;
470: mine->next = PETSC_NULL;
471: if (!packer->rank) packer->n += n;
473: /* add to end of list */
474: if (!next) {
475: packer->next = mine;
476: } else {
477: while (next->next) next = next->next;
478: next->next = mine;
479: }
480: packer->nredundant++;
481: return(0);
482: }
484: /*@C
485: VecPackAddDA - adds a DA vector to a VecPack
487: Collective on VecPack
489: Input Parameter:
490: + packer - the packer object
491: - da - the DA object
492:
493: Level: advanced
495: .seealso VecPackDestroy(), VecPackGather(), VecPackAddDA(), VecPackCreateGlobalVector(),
496: VecPackScatter(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()
498: @*/
499: int VecPackAddDA(VecPack packer,DA da)
500: {
501: int ierr,n;
502: struct VecPackLink *mine,*next = packer->next;
503: Vec global;
506: if (packer->globalvector) {
507: SETERRQ(1,"Cannot add a DA once you have called VecPackCreateGlobalVector()");
508: }
510: /* create new link */
511: PetscNew(struct VecPackLink,&mine);
512: PetscObjectReference((PetscObject)da);
513: DAGetGlobalVector(da,&global);
514: VecGetLocalSize(global,&n);
515: DARestoreGlobalVector(da,&global);
516: mine->n = n;
517: mine->da = da;
518: mine->type = VECPACK_DA;
519: mine->next = PETSC_NULL;
520: packer->n += n;
522: /* add to end of list */
523: if (!next) {
524: packer->next = mine;
525: } else {
526: while (next->next) next = next->next;
527: next->next = mine;
528: }
529: packer->nDA++;
530: return(0);
531: }
533: /*@C
534: VecPackCreateGlobalVector - Creates a vector of the correct size to be gathered into
535: by the packer.
537: Collective on VecPack
539: Input Parameter:
540: . packer - the packer object
542: Output Parameters:
543: . gvec - the global vector
545: Level: advanced
547: Notes: Once this has been created you cannot add additional arrays or vectors to be packed.
549: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackScatter(),
550: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()
552: @*/
553: int VecPackCreateGlobalVector(VecPack packer,Vec *gvec)
554: {
555: int ierr,nprev = 0,rank;
556: struct VecPackLink *next = packer->next;
559: if (packer->globalvector) {
560: VecDuplicate(packer->globalvector,gvec);
561: } else {
562: VecCreateMPI(packer->comm,packer->n,PETSC_DETERMINE,gvec);
563: PetscObjectReference((PetscObject)*gvec);
564: packer->globalvector = *gvec;
566: VecGetSize(*gvec,&packer->N);
567: VecGetOwnershipRange(*gvec,&packer->rstart,PETSC_NULL);
568:
569: /* now set the rstart for each linked array/vector */
570: MPI_Comm_rank(packer->comm,&rank);
571: while (next) {
572: next->rstart = nprev;
573: if (!rank || next->type != VECPACK_ARRAY) nprev += next->n;
574: next = next->next;
575: }
576: }
577: return(0);
578: }
580: /*@C
581: VecPackGetGlobalIndices - Gets the global indices for all the entries in the packed
582: vectors.
584: Collective on VecPack
586: Input Parameter:
587: . packer - the packer object
589: Output Parameters:
590: . idx - the individual indices for each packed vector/array
591:
592: Level: advanced
594: Notes:
595: The idx parameters should be freed by the calling routine with PetscFree()
597: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
598: VecPackGather(), VecPackCreate(), VecPackGetAccess()
600: @*/
601: int VecPackGetGlobalIndices(VecPack packer,...)
602: {
603: va_list Argp;
604: int ierr,i,**idx,n;
605: struct VecPackLink *next = packer->next;
606: Vec global,dglobal;
607: PF pf;
608: PetscScalar *array;
611: VecPackCreateGlobalVector(packer,&global);
613: /* put 0 to N-1 into the global vector */
614: PFCreate(PETSC_COMM_WORLD,1,1,&pf);
615: PFSetType(pf,PFIDENTITY,PETSC_NULL);
616: PFApplyVec(pf,PETSC_NULL,global);
617: PFDestroy(pf);
619: /* loop over packed objects, handling one at at time */
620: va_start(Argp,packer);
621: while (next) {
622: idx = va_arg(Argp, int**);
624: if (next->type == VECPACK_ARRAY) {
625:
626: PetscMalloc(next->n*sizeof(int),idx);
627: if (!packer->rank) {
628: ierr = VecGetArray(global,&array);
629: array += next->rstart;
630: for (i=0; i<next->n; i++) (*idx)[i] = (int)PetscRealPart(array[i]);
631: array -= next->rstart;
632: ierr = VecRestoreArray(global,&array);
633: }
634: MPI_Bcast(*idx,next->n,MPI_INT,0,packer->comm);
636: } else if (next->type == VECPACK_DA) {
637: Vec local;
639: ierr = DACreateLocalVector(next->da,&local);
640: ierr = VecGetArray(global,&array);
641: array += next->rstart;
642: ierr = DAGetGlobalVector(next->da,&dglobal);
643: ierr = VecPlaceArray(dglobal,array);
644: ierr = DAGlobalToLocalBegin(next->da,dglobal,INSERT_VALUES,local);
645: ierr = DAGlobalToLocalEnd(next->da,dglobal,INSERT_VALUES,local);
646: array -= next->rstart;
647: ierr = VecRestoreArray(global,&array);
648: ierr = VecResetArray(dglobal);
649: ierr = DARestoreGlobalVector(next->da,&dglobal);
651: ierr = VecGetArray(local,&array);
652: ierr = VecGetSize(local,&n);
653: ierr = PetscMalloc(n*sizeof(int),idx);
654: for (i=0; i<n; i++) (*idx)[i] = (int)PetscRealPart(array[i]);
655: ierr = VecRestoreArray(local,&array);
656: ierr = VecDestroy(local);
658: } else {
659: SETERRQ(1,"Cannot handle that object type yet");
660: }
661: next = next->next;
662: }
663: va_end(Argp);
664: VecDestroy(global);
665: return(0);
666: }
668: /* -------------------------------------------------------------------------------------*/
669: int VecPackGetLocalVectors_Array(VecPack packer,struct VecPackLink *mine,PetscScalar **array)
670: {
674: PetscMalloc(mine->n*sizeof(PetscScalar),array);
675: return(0);
676: }
678: int VecPackGetLocalVectors_DA(VecPack packer,struct VecPackLink *mine,Vec *local)
679: {
680: int ierr;
682: DAGetLocalVector(mine->da,local);
683: return(0);
684: }
686: int VecPackRestoreLocalVectors_Array(VecPack packer,struct VecPackLink *mine,PetscScalar **array)
687: {
690: PetscFree(*array);
691: return(0);
692: }
694: int VecPackRestoreLocalVectors_DA(VecPack packer,struct VecPackLink *mine,Vec *local)
695: {
696: int ierr;
698: DARestoreLocalVector(mine->da,local);
699: return(0);
700: }
702: /*@C
703: VecPackGetLocalVectors - Gets local vectors and arrays for each part of a VecPack.'
704: Use VecPakcRestoreLocalVectors() to return them.
706: Collective on VecPack
708: Input Parameter:
709: . packer - the packer object
710:
711: Output Parameter:
712: . ... - the individual sequential objects (arrays or vectors)
713:
714: Level: advanced
716: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
717: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
718: VecPackRestoreLocalVectors()
720: @*/
721: int VecPackGetLocalVectors(VecPack packer,...)
722: {
723: va_list Argp;
724: int ierr;
725: struct VecPackLink *next = packer->next;
729: /* loop over packed objects, handling one at at time */
730: va_start(Argp,packer);
731: while (next) {
732: if (next->type == VECPACK_ARRAY) {
733: PetscScalar **array;
734: array = va_arg(Argp, PetscScalar**);
735: VecPackGetLocalVectors_Array(packer,next,array);
736: } else if (next->type == VECPACK_DA) {
737: Vec *vec;
738: vec = va_arg(Argp, Vec*);
739: VecPackGetLocalVectors_DA(packer,next,vec);
740: } else {
741: SETERRQ(1,"Cannot handle that object type yet");
742: }
743: next = next->next;
744: }
745: va_end(Argp);
746: return(0);
747: }
749: /*@C
750: VecPackRestoreLocalVectors - Restores local vectors and arrays for each part of a VecPack.'
751: Use VecPakcRestoreLocalVectors() to return them.
753: Collective on VecPack
755: Input Parameter:
756: . packer - the packer object
757:
758: Output Parameter:
759: . ... - the individual sequential objects (arrays or vectors)
760:
761: Level: advanced
763: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
764: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
765: VecPackGetLocalVectors()
767: @*/
768: int VecPackRestoreLocalVectors(VecPack packer,...)
769: {
770: va_list Argp;
771: int ierr;
772: struct VecPackLink *next = packer->next;
776: /* loop over packed objects, handling one at at time */
777: va_start(Argp,packer);
778: while (next) {
779: if (next->type == VECPACK_ARRAY) {
780: PetscScalar **array;
781: array = va_arg(Argp, PetscScalar**);
782: VecPackRestoreLocalVectors_Array(packer,next,array);
783: } else if (next->type == VECPACK_DA) {
784: Vec *vec;
785: vec = va_arg(Argp, Vec*);
786: VecPackRestoreLocalVectors_DA(packer,next,vec);
787: } else {
788: SETERRQ(1,"Cannot handle that object type yet");
789: }
790: next = next->next;
791: }
792: va_end(Argp);
793: return(0);
794: }
796: /* -------------------------------------------------------------------------------------*/
797: int VecPackGetEntries_Array(VecPack packer,struct VecPackLink *mine,int *n)
798: {
800: if (n) *n = mine->n;
801: return(0);
802: }
804: int VecPackGetEntries_DA(VecPack packer,struct VecPackLink *mine,DA *da)
805: {
807: if (da) *da = mine->da;
808: return(0);
809: }
811: /*@C
812: VecPackGetEntries - Gets the DA, redundant size, etc for each entry in a VecPack.
813: Use VecPackRestoreEntries() to return them.
815: Collective on VecPack
817: Input Parameter:
818: . packer - the packer object
819:
820: Output Parameter:
821: . ... - the individual entries, DAs or integer sizes)
822:
823: Level: advanced
825: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
826: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(),
827: VecPackRestoreLocalVectors(), VecPackGetLocalVectors(), VecPackRestoreEntries()
829: @*/
830: int VecPackGetEntries(VecPack packer,...)
831: {
832: va_list Argp;
833: int ierr;
834: struct VecPackLink *next = packer->next;
838: /* loop over packed objects, handling one at at time */
839: va_start(Argp,packer);
840: while (next) {
841: if (next->type == VECPACK_ARRAY) {
842: int *n;
843: n = va_arg(Argp, int*);
844: VecPackGetEntries_Array(packer,next,n);
845: } else if (next->type == VECPACK_DA) {
846: DA *da;
847: da = va_arg(Argp, DA*);
848: VecPackGetEntries_DA(packer,next,da);
849: } else {
850: SETERRQ(1,"Cannot handle that object type yet");
851: }
852: next = next->next;
853: }
854: va_end(Argp);
855: return(0);
856: }
858: /*@C
859: VecPackRefine - Refines a VecPack by refining all of its DAs
861: Collective on VecPack
863: Input Parameters:
864: + packer - the packer object
865: - comm - communicator to contain the new DM object, usually PETSC_NULL
867: Output Parameter:
868: . fine - new packer
869:
870: Level: advanced
872: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
873: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()
875: @*/
876: int VecPackRefine(VecPack packer,MPI_Comm comm,VecPack *fine)
877: {
878: int ierr;
879: struct VecPackLink *next = packer->next;
880: DA da;
883: VecPackCreate(comm,fine);
885: /* loop over packed objects, handling one at at time */
886: while (next) {
887: if (next->type == VECPACK_ARRAY) {
888: VecPackAddArray(*fine,next->n);
889: } else if (next->type == VECPACK_DA) {
890: DARefine(next->da,comm,&da);
891: VecPackAddDA(*fine,da);
892: PetscObjectDereference((PetscObject)da);
893: } else {
894: SETERRQ(1,"Cannot handle that object type yet");
895: }
896: next = next->next;
897: }
898: return(0);
899: }
901: #include petscmat.h
903: struct MatPackLink {
904: Mat A;
905: struct MatPackLink *next;
906: };
908: struct MatPack {
909: VecPack right,left;
910: struct MatPackLink *next;
911: };
913: int MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscTruth add)
914: {
915: struct MatPack *mpack;
916: struct VecPackLink *xnext,*ynext;
917: struct MatPackLink *anext;
918: PetscScalar *xarray,*yarray;
919: int ierr,i;
920: Vec xglobal,yglobal;
923: MatShellGetContext(A,(void**)&mpack);
924: xnext = mpack->right->next;
925: ynext = mpack->left->next;
926: anext = mpack->next;
928: while (xnext) {
929: if (xnext->type == VECPACK_ARRAY) {
930: if (!mpack->right->rank) {
931: ierr = VecGetArray(x,&xarray);
932: ierr = VecGetArray(y,&yarray);
933: if (add) {
934: for (i=0; i<xnext->n; i++) {
935: yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
936: }
937: } else {
938: ierr = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));
939: }
940: ierr = VecRestoreArray(x,&xarray);
941: ierr = VecRestoreArray(y,&yarray);
942: }
943: } else if (xnext->type == VECPACK_DA) {
944: ierr = VecGetArray(x,&xarray);
945: ierr = VecGetArray(y,&yarray);
946: ierr = DAGetGlobalVector(xnext->da,&xglobal);
947: ierr = DAGetGlobalVector(ynext->da,&yglobal);
948: ierr = VecPlaceArray(xglobal,xarray+xnext->rstart);
949: ierr = VecPlaceArray(yglobal,yarray+ynext->rstart);
950: if (add) {
951: ierr = MatMultAdd(anext->A,xglobal,yglobal,yglobal);
952: } else {
953: ierr = MatMult(anext->A,xglobal,yglobal);
954: }
955: ierr = VecRestoreArray(x,&xarray);
956: ierr = VecRestoreArray(y,&yarray);
957: ierr = VecResetArray(xglobal);
958: ierr = VecResetArray(yglobal);
959: ierr = DARestoreGlobalVector(xnext->da,&xglobal);
960: ierr = DARestoreGlobalVector(ynext->da,&yglobal);
961: anext = anext->next;
962: } else {
963: SETERRQ(1,"Cannot handle that object type yet");
964: }
965: xnext = xnext->next;
966: ynext = ynext->next;
967: }
968: return(0);
969: }
971: int MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
972: {
975: if (z != y) SETERRQ(1,"Handles y == z only");
976: MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);
977: return(0);
978: }
980: int MatMult_Shell_Pack(Mat A,Vec x,Vec y)
981: {
984: MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);
985: return(0);
986: }
988: int MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
989: {
990: struct MatPack *mpack;
991: struct VecPackLink *xnext,*ynext;
992: struct MatPackLink *anext;
993: PetscScalar *xarray,*yarray;
994: int ierr;
995: Vec xglobal,yglobal;
998: ierr = MatShellGetContext(A,(void**)&mpack);
999: xnext = mpack->left->next;
1000: ynext = mpack->right->next;
1001: anext = mpack->next;
1003: while (xnext) {
1004: if (xnext->type == VECPACK_ARRAY) {
1005: if (!mpack->right->rank) {
1006: ierr = VecGetArray(x,&xarray);
1007: ierr = VecGetArray(y,&yarray);
1008: ierr = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));
1009: ierr = VecRestoreArray(x,&xarray);
1010: ierr = VecRestoreArray(y,&yarray);
1011: }
1012: } else if (xnext->type == VECPACK_DA) {
1013: ierr = VecGetArray(x,&xarray);
1014: ierr = VecGetArray(y,&yarray);
1015: ierr = DAGetGlobalVector(xnext->da,&xglobal);
1016: ierr = DAGetGlobalVector(ynext->da,&yglobal);
1017: ierr = VecPlaceArray(xglobal,xarray+xnext->rstart);
1018: ierr = VecPlaceArray(yglobal,yarray+ynext->rstart);
1019: ierr = MatMultTranspose(anext->A,xglobal,yglobal);
1020: ierr = VecRestoreArray(x,&xarray);
1021: ierr = VecRestoreArray(y,&yarray);
1022: ierr = VecResetArray(xglobal);
1023: ierr = VecResetArray(yglobal);
1024: ierr = DARestoreGlobalVector(xnext->da,&xglobal);
1025: ierr = DARestoreGlobalVector(ynext->da,&yglobal);
1026: anext = anext->next;
1027: } else {
1028: SETERRQ(1,"Cannot handle that object type yet");
1029: }
1030: xnext = xnext->next;
1031: ynext = ynext->next;
1032: }
1033: return(0);
1034: }
1036: int MatDestroy_Shell_Pack(Mat A)
1037: {
1038: struct MatPack *mpack;
1039: struct MatPackLink *anext,*oldanext;
1040: int ierr;
1043: ierr = MatShellGetContext(A,(void**)&mpack);
1044: anext = mpack->next;
1046: while (anext) {
1047: ierr = MatDestroy(anext->A);
1048: oldanext = anext;
1049: anext = anext->next;
1050: ierr = PetscFree(oldanext);
1051: }
1052: PetscFree(mpack);
1053: return(0);
1054: }
1056: /*@C
1057: VecPackGetInterpolation - GetInterpolations a VecPack by refining all of its DAs
1059: Collective on VecPack
1061: Input Parameters:
1062: + coarse - coarse grid packer
1063: - fine - fine grid packer
1065: Output Parameter:
1066: + A - interpolation matrix
1067: - v - scaling vector
1068:
1069: Level: advanced
1071: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
1072: VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()
1074: @*/
1075: int VecPackGetInterpolation(VecPack coarse,VecPack fine,Mat *A,Vec *v)
1076: {
1077: int ierr,m,n,M,N;
1078: struct VecPackLink *nextc = coarse->next;
1079: struct VecPackLink *nextf = fine->next;
1080: struct MatPackLink *nextmat,*pnextmat = 0;
1081: struct MatPack *mpack;
1082: Vec gcoarse,gfine;
1085: /* use global vectors only for determining matrix layout */
1086: VecPackCreateGlobalVector(coarse,&gcoarse);
1087: VecPackCreateGlobalVector(fine,&gfine);
1088: VecGetLocalSize(gcoarse,&n);
1089: VecGetLocalSize(gfine,&m);
1090: VecGetSize(gcoarse,&N);
1091: VecGetSize(gfine,&M);
1092: VecDestroy(gcoarse);
1093: VecDestroy(gfine);
1095: ierr = PetscNew(struct MatPack,&mpack);
1096: mpack->right = coarse;
1097: mpack->left = fine;
1098: ierr = MatCreateShell(fine->comm,m,n,M,N,mpack,A);
1099: ierr = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);
1100: ierr = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);
1101: ierr = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);
1102: ierr = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);
1104: /* loop over packed objects, handling one at at time */
1105: while (nextc) {
1106: if (nextc->type != nextf->type) SETERRQ(1,"Two VecPack have different layout");
1108: if (nextc->type == VECPACK_ARRAY) {
1109: ;
1110: } else if (nextc->type == VECPACK_DA) {
1111: ierr = PetscNew(struct MatPackLink,&nextmat);
1112: nextmat->next = 0;
1113: if (pnextmat) {
1114: pnextmat->next = nextmat;
1115: pnextmat = nextmat;
1116: } else {
1117: pnextmat = nextmat;
1118: mpack->next = nextmat;
1119: }
1120: DAGetInterpolation(nextc->da,nextf->da,&nextmat->A,PETSC_NULL);
1121: } else {
1122: SETERRQ(1,"Cannot handle that object type yet");
1123: }
1124: nextc = nextc->next;
1125: nextf = nextf->next;
1126: }
1127: return(0);
1128: }