Actual source code: pack.c

  1: /*$Id: pack.c,v 1.16 2001/03/28 19:42:48 balay Exp $*/
  2: 
  3: #include "petscda.h"     /*I      "petscda.h"     I*/
  4: #include "petscmat.h"    /*I      "petscmat.h"    I*/

  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,ISColoring*,Mat*);
 25:   int  (*getinterpolation)(VecPack,VecPack,Mat*,Vec*);
 26:   int  (*refine)(VecPack,MPI_Comm,VecPack*);
 27: };

 29: struct _p_VecPack {
 30:   PETSCHEADER(struct _VecPackOps)
 31:   int                rank;
 32:   int                n,N,rstart;   /* rstart is relative to all processors */
 33:   Vec                globalvector;
 34:   int                nDA,nredundant;
 35:   struct VecPackLink *next;
 36: };

 38: /*@C
 39:     VecPackCreate - Creates a vector packer, used to generate "composite"
 40:       vectors made up of several subvectors.

 42:     Collective on MPI_Comm

 44:     Input Parameter:
 45: .   comm - the processors that will share the global vector

 47:     Output Parameters:
 48: .   packer - the packer object

 50:     Level: advanced

 52: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackScatter(),
 53:          VecPackGather(), VecPackCreateGlobalVector(), VecPackGetGlobalIndices(), VecPackGetAccess()

 55: @*/
 56: int VecPackCreate(MPI_Comm comm,VecPack *packer)
 57: {
 58:   int     ierr;
 59:   VecPack p;

 62:   PetscHeaderCreate(p,_p_VecPack,struct _VecPackOps,DA_COOKIE,0,"VecPack",comm,VecPackDestroy,0);
 63:   PetscLogObjectCreate(p);
 64:   p->n            = 0;
 65:   p->next         = PETSC_NULL;
 66:   p->comm         = comm;
 67:   p->globalvector = PETSC_NULL;
 68:   p->nredundant   = 0;
 69:   p->nDA          = 0;
 70:   ierr            = MPI_Comm_rank(comm,&p->rank);

 72:   p->ops->createglobalvector = VecPackCreateGlobalVector;
 73:   p->ops->refine             = VecPackRefine;
 74:   p->ops->getinterpolation   = VecPackGetInterpolation;
 75:   *packer = p;
 76:   return(0);
 77: }

 79: /*@C
 80:     VecPackDestroy - Destroys a vector packer.

 82:     Collective on VecPack

 84:     Input Parameter:
 85: .   packer - the packer object

 87:     Level: advanced

 89: .seealso VecPackCreate(), VecPackAddArray(), VecPackAddDA(), VecPackScatter(),
 90:          VecPackGather(), VecPackCreateGlobalVector(), VecPackGetGlobalIndices(), VecPackGetAccess()

 92: @*/
 93: int VecPackDestroy(VecPack packer)
 94: {
 95:   int                ierr;
 96:   struct VecPackLink *next = packer->next,*prev;

 99:   if (--packer->refct > 0) return(0);
100:   while (next) {
101:     prev = next;
102:     next = next->next;
103:     if (prev->type == VECPACK_DA) {
104:       DADestroy(prev->da);
105:     }
106:     PetscFree(prev);
107:   }
108:   if (packer->globalvector) {
109:     VecDestroy(packer->globalvector);
110:   }
111:   PetscHeaderDestroy(packer);
112:   return(0);
113: }

115: /* --------------------------------------------------------------------------------------*/

117: int VecPackGetAccess_Array(VecPack packer,struct VecPackLink *mine,Vec vec,Scalar **array)
118: {
119:   int    ierr;
120:   Scalar *varray;

123:   if (array) {
124:     if (!packer->rank) {
125:       ierr    = VecGetArray(vec,&varray);
126:       *array  = varray + mine->rstart;
127:       ierr    = VecRestoreArray(vec,&varray);
128:     } else {
129:       *array = 0;
130:     }
131:   }
132:   return(0);
133: }

135: int VecPackGetAccess_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec *global)
136: {
137:   int    ierr;
138:   Scalar *array;

141:   if (global) {
142:     ierr    = DAGetGlobalVector(mine->da,global);
143:     ierr    = VecGetArray(vec,&array);
144:     ierr    = VecPlaceArray(*global,array+mine->rstart);
145:     ierr    = VecRestoreArray(vec,&array);
146:   }
147:   return(0);
148: }

150: int VecPackRestoreAccess_Array(VecPack packer,struct VecPackLink *mine,Vec vec,Scalar **array)
151: {
153:   return(0);
154: }

156: int VecPackRestoreAccess_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec *global)
157: {
158:   int    ierr;

161:   if (global) {
162:     VecResetArray(*global);
163:     DARestoreGlobalVector(mine->da,global);
164:   }
165:   return(0);
166: }

168: int VecPackScatter_Array(VecPack packer,struct VecPackLink *mine,Vec vec,Scalar *array)
169: {
170:   int    ierr;
171:   Scalar *varray;


175:   if (!packer->rank) {
176:     ierr    = VecGetArray(vec,&varray);
177:     ierr    = PetscMemcpy(array,varray+mine->rstart,mine->n*sizeof(Scalar));
178:     ierr    = VecRestoreArray(vec,&varray);
179:   }
180:   ierr    = MPI_Bcast(array,mine->n,MPIU_SCALAR,0,packer->comm);
181:   return(0);
182: }

184: int VecPackScatter_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec local)
185: {
186:   int    ierr;
187:   Scalar *array;
188:   Vec    global;

191:   DAGetGlobalVector(mine->da,&global);
192:   VecGetArray(vec,&array);
193:   VecPlaceArray(global,array+mine->rstart);
194:   DAGlobalToLocalBegin(mine->da,global,INSERT_VALUES,local);
195:   DAGlobalToLocalEnd(mine->da,global,INSERT_VALUES,local);
196:   VecRestoreArray(vec,&array);
197:   VecResetArray(global);
198:   DARestoreGlobalVector(mine->da,&global);
199:   return(0);
200: }

202: int VecPackGather_Array(VecPack packer,struct VecPackLink *mine,Vec vec,Scalar *array)
203: {
204:   int    ierr;
205:   Scalar *varray;

208:   if (!packer->rank) {
209:     ierr    = VecGetArray(vec,&varray);
210:     if (varray+mine->rstart == array) SETERRQ(1,"You need not VecPackGather() into objects obtained via VecPackGetAccess()");
211:     ierr    = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(Scalar));
212:     ierr    = VecRestoreArray(vec,&varray);
213:   }
214:   return(0);
215: }

217: int VecPackGather_DA(VecPack packer,struct VecPackLink *mine,Vec vec,Vec local)
218: {
219:   int    ierr;
220:   Scalar *array;
221:   Vec    global;

224:   DAGetGlobalVector(mine->da,&global);
225:   VecGetArray(vec,&array);
226:   VecPlaceArray(global,array+mine->rstart);
227:   DALocalToGlobal(mine->da,local,INSERT_VALUES,global);
228:   VecRestoreArray(vec,&array);
229:   VecResetArray(global);
230:   DARestoreGlobalVector(mine->da,&global);
231:   return(0);
232: }

234: /* ----------------------------------------------------------------------------------*/

236: #include <stdarg.h>

238: /*@C
239:     VecPackGetAccess - Allows one to access the individual packed vectors in their global
240:        representation.

242:     Collective on VecPack

244:     Input Parameter:
245: +    packer - the packer object
246: .    gvec - the global vector
247: -    ... - the individual sequential or parallel objects (arrays or vectors)
248:  
249:     Level: advanced

251: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
252:          VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackScatter(),
253:          VecPackRestoreAccess()

255: @*/
256: int VecPackGetAccess(VecPack packer,Vec gvec,...)
257: {
258:   va_list            Argp;
259:   int                ierr;
260:   struct VecPackLink *next = packer->next;

263:   if (!packer->globalvector) {
264:     SETERRQ(1,"Must first create global vector with VecPackCreateGlobalVector()");
265:   }

267:   /* loop over packed objects, handling one at at time */
268:   va_start(Argp,gvec);
269:   while (next) {
270:     if (next->type == VECPACK_ARRAY) {
271:       Scalar **array;
272:       array = va_arg(Argp, Scalar**);
273:       ierr  = VecPackGetAccess_Array(packer,next,gvec,array);
274:     } else if (next->type == VECPACK_DA) {
275:       Vec *vec;
276:       vec  = va_arg(Argp, Vec*);
277:       VecPackGetAccess_DA(packer,next,gvec,vec);
278:     } else {
279:       SETERRQ(1,"Cannot handle that object type yet");
280:     }
281:     next = next->next;
282:   }
283:   va_end(Argp);
284:   return(0);
285: }

287: /*@C
288:     VecPackRestoreAccess - Allows one to access the individual packed vectors in their global
289:        representation.

291:     Collective on VecPack

293:     Input Parameter:
294: +    packer - the packer object
295: .    gvec - the global vector
296: -    ... - the individual sequential or parallel objects (arrays or vectors)
297:  
298:     Level: advanced

300: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
301:          VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackScatter(),
302:          VecPackRestoreAccess()

304: @*/
305: int VecPackRestoreAccess(VecPack packer,Vec gvec,...)
306: {
307:   va_list            Argp;
308:   int                ierr;
309:   struct VecPackLink *next = packer->next;

312:   if (!packer->globalvector) {
313:     SETERRQ(1,"Must first create global vector with VecPackCreateGlobalVector()");
314:   }

316:   /* loop over packed objects, handling one at at time */
317:   va_start(Argp,gvec);
318:   while (next) {
319:     if (next->type == VECPACK_ARRAY) {
320:       Scalar **array;
321:       array = va_arg(Argp, Scalar**);
322:       ierr  = VecPackRestoreAccess_Array(packer,next,gvec,array);
323:     } else if (next->type == VECPACK_DA) {
324:       Vec *vec;
325:       vec  = va_arg(Argp, Vec*);
326:       VecPackRestoreAccess_DA(packer,next,gvec,vec);
327:     } else {
328:       SETERRQ(1,"Cannot handle that object type yet");
329:     }
330:     next = next->next;
331:   }
332:   va_end(Argp);
333:   return(0);
334: }

336: /*@C
337:     VecPackScatter - Scatters from a global packed vector into its individual local vectors

339:     Collective on VecPack

341:     Input Parameter:
342: +    packer - the packer object
343: .    gvec - the global vector
344: -    ... - the individual sequential objects (arrays or vectors)
345:  
346:     Level: advanced

348: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
349:          VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()

351: @*/
352: int VecPackScatter(VecPack packer,Vec gvec,...)
353: {
354:   va_list            Argp;
355:   int                ierr;
356:   struct VecPackLink *next = packer->next;

359:   if (!packer->globalvector) {
360:     SETERRQ(1,"Must first create global vector with VecPackCreateGlobalVector()");
361:   }

363:   /* loop over packed objects, handling one at at time */
364:   va_start(Argp,gvec);
365:   while (next) {
366:     if (next->type == VECPACK_ARRAY) {
367:       Scalar *array;
368:       array = va_arg(Argp, Scalar*);
369:       VecPackScatter_Array(packer,next,gvec,array);
370:     } else if (next->type == VECPACK_DA) {
371:       Vec vec;
372:       vec = va_arg(Argp, Vec);
374:       VecPackScatter_DA(packer,next,gvec,vec);
375:     } else {
376:       SETERRQ(1,"Cannot handle that object type yet");
377:     }
378:     next = next->next;
379:   }
380:   va_end(Argp);
381:   return(0);
382: }

384: /*@C
385:     VecPackGather - Gathers into a global packed vector from its individual local vectors

387:     Collective on VecPack

389:     Input Parameter:
390: +    packer - the packer object
391: .    gvec - the global vector
392: -    ... - the individual sequential objects (arrays or vectors)
393:  
394:     Level: advanced

396: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
397:          VecPackScatter(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()

399: @*/
400: int VecPackGather(VecPack packer,Vec gvec,...)
401: {
402:   va_list            Argp;
403:   int                ierr;
404:   struct VecPackLink *next = packer->next;

407:   if (!packer->globalvector) {
408:     SETERRQ(1,"Must first create global vector with VecPackCreateGlobalVector()");
409:   }

411:   /* loop over packed objects, handling one at at time */
412:   va_start(Argp,gvec);
413:   while (next) {
414:     if (next->type == VECPACK_ARRAY) {
415:       Scalar *array;
416:       array = va_arg(Argp, Scalar*);
417:       ierr  = VecPackGather_Array(packer,next,gvec,array);
418:     } else if (next->type == VECPACK_DA) {
419:       Vec vec;
420:       vec = va_arg(Argp, Vec);
422:       VecPackGather_DA(packer,next,gvec,vec);
423:     } else {
424:       SETERRQ(1,"Cannot handle that object type yet");
425:     }
426:     next = next->next;
427:   }
428:   va_end(Argp);
429:   return(0);
430: }

432: /*@C
433:     VecPackAddArray - adds an "redundant" array to a VecPack. The array values will 
434:        be stored in part of the array on processor 0.

436:     Collective on VecPack

438:     Input Parameter:
439: +    packer - the packer object
440: -    n - the length of the array
441:  
442:     Level: advanced

444: .seealso VecPackDestroy(), VecPackGather(), VecPackAddDA(), VecPackCreateGlobalVector(),
445:          VecPackScatter(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()

447: @*/
448: int VecPackAddArray(VecPack packer,int n)
449: {
450:   struct VecPackLink *mine,*next = packer->next;

454:   if (packer->globalvector) {
455:     SETERRQ(1,"Cannot add an array once you have called VecPackCreateGlobalVector()");
456:   }

458:   /* create new link */
459:   ierr                = PetscNew(struct VecPackLink,&mine);
460:   mine->n             = n;
461:   mine->da            = PETSC_NULL;
462:   mine->type          = VECPACK_ARRAY;
463:   mine->next          = PETSC_NULL;
464:   if (!packer->rank) packer->n += n;

466:   /* add to end of list */
467:   if (!next) {
468:     packer->next = mine;
469:   } else {
470:     while (next->next) next = next->next;
471:     next->next = mine;
472:   }
473:   packer->nredundant++;
474:   return(0);
475: }

477: /*@C
478:     VecPackAddDA - adds a DA vector to a VecPack

480:     Collective on VecPack

482:     Input Parameter:
483: +    packer - the packer object
484: -    da - the DA object
485:  
486:     Level: advanced

488: .seealso VecPackDestroy(), VecPackGather(), VecPackAddDA(), VecPackCreateGlobalVector(),
489:          VecPackScatter(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()

491: @*/
492: int VecPackAddDA(VecPack packer,DA da)
493: {
494:   int                ierr,n;
495:   struct VecPackLink *mine,*next = packer->next;
496:   Vec                global;

499:   if (packer->globalvector) {
500:     SETERRQ(1,"Cannot add a DA once you have called VecPackCreateGlobalVector()");
501:   }

503:   /* create new link */
504:   PetscNew(struct VecPackLink,&mine);
505:   PetscObjectReference((PetscObject)da);
506:   DAGetGlobalVector(da,&global);
507:   VecGetLocalSize(global,&n);
508:   DARestoreGlobalVector(da,&global);
509:   mine->n      = n;
510:   mine->da     = da;
511:   mine->type   = VECPACK_DA;
512:   mine->next   = PETSC_NULL;
513:   packer->n   += n;

515:   /* add to end of list */
516:   if (!next) {
517:     packer->next = mine;
518:   } else {
519:     while (next->next) next = next->next;
520:     next->next = mine;
521:   }
522:   packer->nDA++;
523:   return(0);
524: }

526: /*@C
527:     VecPackCreateGlobalVector - Creates a vector of the correct size to be gathered into 
528:         by the packer.

530:     Collective on VecPack

532:     Input Parameter:
533: .    packer - the packer object

535:     Output Parameters:
536: .   gvec - the global vector

538:     Level: advanced

540:     Notes: Once this has been created you cannot add additional arrays or vectors to be packed.

542: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackScatter(),
543:          VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()

545: @*/
546: int VecPackCreateGlobalVector(VecPack packer,Vec *gvec)
547: {
548:   int                ierr,nprev = 0,rank;
549:   struct VecPackLink *next = packer->next;

552:   if (packer->globalvector) {
553:     VecDuplicate(packer->globalvector,gvec);
554:   } else {
555:     VecCreateMPI(packer->comm,packer->n,PETSC_DETERMINE,gvec);
556:     PetscObjectReference((PetscObject)*gvec);
557:     packer->globalvector = *gvec;

559:     VecGetSize(*gvec,&packer->N);
560:     VecGetOwnershipRange(*gvec,&packer->rstart,PETSC_NULL);
561: 
562:     /* now set the rstart for each linked array/vector */
563:     MPI_Comm_rank(packer->comm,&rank);
564:     while (next) {
565:       next->rstart = nprev;
566:       if (!rank || next->type != VECPACK_ARRAY) nprev += next->n;
567:       next = next->next;
568:     }
569:   }
570:   return(0);
571: }

573: /*@C
574:     VecPackGetGlobalIndices - Gets the global indices for all the entries in the packed
575:       vectors.

577:     Collective on VecPack

579:     Input Parameter:
580: .    packer - the packer object

582:     Output Parameters:
583: .    idx - the individual indices for each packed vector/array
584:  
585:     Level: advanced

587:     Notes:
588:        The idx parameters should be freed by the calling routine with PetscFree()

590: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
591:          VecPackGather(), VecPackCreate(), VecPackGetAccess()

593: @*/
594: int VecPackGetGlobalIndices(VecPack packer,...)
595: {
596:   va_list            Argp;
597:   int                ierr,i,**idx,n;
598:   struct VecPackLink *next = packer->next;
599:   Vec                global,dglobal;
600:   PF                 pf;
601:   Scalar             *array;

604:   VecPackCreateGlobalVector(packer,&global);

606:   /* put 0 to N-1 into the global vector */
607:   PFCreate(PETSC_COMM_WORLD,1,1,&pf);
608:   PFSetType(pf,PFIDENTITY,PETSC_NULL);
609:   PFApplyVec(pf,PETSC_NULL,global);
610:   PFDestroy(pf);

612:   /* loop over packed objects, handling one at at time */
613:   va_start(Argp,packer);
614:   while (next) {
615:     idx = va_arg(Argp, int**);

617:     if (next->type == VECPACK_ARRAY) {
618: 
619:       PetscMalloc(next->n*sizeof(int),idx);
620:       if (!packer->rank) {
621:         ierr   = VecGetArray(global,&array);
622:         array += next->rstart;
623:         for (i=0; i<next->n; i++) (*idx)[i] = (int)PetscRealPart(array[i]);
624:         array -= next->rstart;
625:         ierr   = VecRestoreArray(global,&array);
626:       }
627:       MPI_Bcast(*idx,next->n,MPI_INT,0,packer->comm);

629:     } else if (next->type == VECPACK_DA) {
630:       Vec local;

632:       ierr   = DACreateLocalVector(next->da,&local);
633:       ierr   = VecGetArray(global,&array);
634:       array += next->rstart;
635:       ierr   = DAGetGlobalVector(next->da,&dglobal);
636:       ierr   = VecPlaceArray(dglobal,array);
637:       ierr   = DAGlobalToLocalBegin(next->da,dglobal,INSERT_VALUES,local);
638:       ierr   = DAGlobalToLocalEnd(next->da,dglobal,INSERT_VALUES,local);
639:       array -= next->rstart;
640:       ierr   = VecRestoreArray(global,&array);
641:       ierr   = VecResetArray(dglobal);
642:       ierr   = DARestoreGlobalVector(next->da,&dglobal);

644:       ierr   = VecGetArray(local,&array);
645:       ierr   = VecGetSize(local,&n);
646:       ierr   = PetscMalloc(n*sizeof(int),idx);
647:       for (i=0; i<n; i++) (*idx)[i] = (int)PetscRealPart(array[i]);
648:       ierr    = VecRestoreArray(local,&array);
649:       ierr    = VecDestroy(local);

651:     } else {
652:       SETERRQ(1,"Cannot handle that object type yet");
653:     }
654:     next = next->next;
655:   }
656:   va_end(Argp);
657:   VecDestroy(global);
658:   return(0);
659: }

661: /* -------------------------------------------------------------------------------------*/
662: int VecPackGetLocalVectors_Array(VecPack packer,struct VecPackLink *mine,Scalar **array)
663: {

667:   PetscMalloc(mine->n*sizeof(Scalar),array);
668:   return(0);
669: }

671: int VecPackGetLocalVectors_DA(VecPack packer,struct VecPackLink *mine,Vec *local)
672: {
673:   int    ierr;
675:   DAGetLocalVector(mine->da,local);
676:   return(0);
677: }

679: int VecPackRestoreLocalVectors_Array(VecPack packer,struct VecPackLink *mine,Scalar **array)
680: {
683:   PetscFree(*array);
684:   return(0);
685: }

687: int VecPackRestoreLocalVectors_DA(VecPack packer,struct VecPackLink *mine,Vec *local)
688: {
689:   int    ierr;
691:   DARestoreLocalVector(mine->da,local);
692:   return(0);
693: }

695: /*@C
696:     VecPackGetLocalVectors - Gets local vectors and arrays for each part of a VecPack.'
697:        Use VecPakcRestoreLocalVectors() to return them.

699:     Collective on VecPack

701:     Input Parameter:
702: .    packer - the packer object
703:  
704:     Output Parameter:
705: .    ... - the individual sequential objects (arrays or vectors)
706:  
707:     Level: advanced

709: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
710:          VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(), 
711:          VecPackRestoreLocalVectors()

713: @*/
714: int VecPackGetLocalVectors(VecPack packer,...)
715: {
716:   va_list            Argp;
717:   int                ierr;
718:   struct VecPackLink *next = packer->next;


722:   /* loop over packed objects, handling one at at time */
723:   va_start(Argp,packer);
724:   while (next) {
725:     if (next->type == VECPACK_ARRAY) {
726:       Scalar **array;
727:       array = va_arg(Argp, Scalar**);
728:       VecPackGetLocalVectors_Array(packer,next,array);
729:     } else if (next->type == VECPACK_DA) {
730:       Vec *vec;
731:       vec = va_arg(Argp, Vec*);
732:       VecPackGetLocalVectors_DA(packer,next,vec);
733:     } else {
734:       SETERRQ(1,"Cannot handle that object type yet");
735:     }
736:     next = next->next;
737:   }
738:   va_end(Argp);
739:   return(0);
740: }

742: /*@C
743:     VecPackRestoreLocalVectors - Restores local vectors and arrays for each part of a VecPack.'
744:        Use VecPakcRestoreLocalVectors() to return them.

746:     Collective on VecPack

748:     Input Parameter:
749: .    packer - the packer object
750:  
751:     Output Parameter:
752: .    ... - the individual sequential objects (arrays or vectors)
753:  
754:     Level: advanced

756: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
757:          VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(), 
758:          VecPackGetLocalVectors()

760: @*/
761: int VecPackRestoreLocalVectors(VecPack packer,...)
762: {
763:   va_list            Argp;
764:   int                ierr;
765:   struct VecPackLink *next = packer->next;


769:   /* loop over packed objects, handling one at at time */
770:   va_start(Argp,packer);
771:   while (next) {
772:     if (next->type == VECPACK_ARRAY) {
773:       Scalar **array;
774:       array = va_arg(Argp, Scalar**);
775:       VecPackRestoreLocalVectors_Array(packer,next,array);
776:     } else if (next->type == VECPACK_DA) {
777:       Vec *vec;
778:       vec = va_arg(Argp, Vec*);
779:       VecPackRestoreLocalVectors_DA(packer,next,vec);
780:     } else {
781:       SETERRQ(1,"Cannot handle that object type yet");
782:     }
783:     next = next->next;
784:   }
785:   va_end(Argp);
786:   return(0);
787: }

789: /* -------------------------------------------------------------------------------------*/
790: int VecPackGetEntries_Array(VecPack packer,struct VecPackLink *mine,int *n)
791: {
793:   if (n) *n = mine->n;
794:   return(0);
795: }

797: int VecPackGetEntries_DA(VecPack packer,struct VecPackLink *mine,DA *da)
798: {
800:   if (da) *da = mine->da;
801:   return(0);
802: }

804: /*@C
805:     VecPackGetEntries - Gets the DA, redundant size, etc for each entry in a VecPack.
806:        Use VecPakcRestoreEntries() to return them.

808:     Collective on VecPack

810:     Input Parameter:
811: .    packer - the packer object
812:  
813:     Output Parameter:
814: .    ... - the individual entries, DAs or integer sizes)
815:  
816:     Level: advanced

818: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
819:          VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess(), 
820:          VecPackRestoreLocalVectors(), VecPackGetLocalVectors()

822: @*/
823: int VecPackGetEntries(VecPack packer,...)
824: {
825:   va_list            Argp;
826:   int                ierr;
827:   struct VecPackLink *next = packer->next;


831:   /* loop over packed objects, handling one at at time */
832:   va_start(Argp,packer);
833:   while (next) {
834:     if (next->type == VECPACK_ARRAY) {
835:       int *n;
836:       n = va_arg(Argp, int*);
837:       VecPackGetEntries_Array(packer,next,n);
838:     } else if (next->type == VECPACK_DA) {
839:       DA *da;
840:       da = va_arg(Argp, DA*);
841:       VecPackGetEntries_DA(packer,next,da);
842:     } else {
843:       SETERRQ(1,"Cannot handle that object type yet");
844:     }
845:     next = next->next;
846:   }
847:   va_end(Argp);
848:   return(0);
849: }

851: /*@C
852:     VecPackRefine - Refines a VecPack by refining all of its DAs

854:     Collective on VecPack

856:     Input Parameters:
857: +    packer - the packer object
858: -    comm - communicator to contain the new DM object, usually PETSC_NULL

860:     Output Parameter:
861: .    fine - new packer
862:  
863:     Level: advanced

865: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
866:          VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()

868: @*/
869: int VecPackRefine(VecPack packer,MPI_Comm comm,VecPack *fine)
870: {
871:   int                ierr;
872:   struct VecPackLink *next = packer->next;
873:   DA                 da;

876:   VecPackCreate(comm,fine);

878:   /* loop over packed objects, handling one at at time */
879:   while (next) {
880:     if (next->type == VECPACK_ARRAY) {
881:       VecPackAddArray(*fine,next->n);
882:     } else if (next->type == VECPACK_DA) {
883:       DARefine(next->da,comm,&da);
884:       VecPackAddDA(*fine,da);
885:       PetscObjectDereference((PetscObject)da);
886:     } else {
887:       SETERRQ(1,"Cannot handle that object type yet");
888:     }
889:     next = next->next;
890:   }
891:   return(0);
892: }

894: #include "petscmat.h"

896: struct MatPackLink {
897:   Mat                A;
898:   struct MatPackLink *next;
899: };

901: struct MatPack {
902:   VecPack            right,left;
903:   struct MatPackLink *next;
904: };

906: int MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscTruth add)
907: {
908:   struct MatPack     *mpack;
909:   struct VecPackLink *xnext,*ynext;
910:   struct MatPackLink *anext;
911:   Scalar             *xarray,*yarray;
912:   int                ierr,i;
913:   Vec                xglobal,yglobal;

916:   MatShellGetContext(A,(void**)&mpack);
917:   xnext = mpack->right->next;
918:   ynext = mpack->left->next;
919:   anext = mpack->next;

921:   while (xnext) {
922:     if (xnext->type == VECPACK_ARRAY) {
923:       if (!mpack->right->rank) {
924:         ierr    = VecGetArray(x,&xarray);
925:         ierr    = VecGetArray(y,&yarray);
926:         if (add) {
927:           for (i=0; i<xnext->n; i++) {
928:             yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
929:           }
930:         } else {
931:           ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(Scalar));
932:         }
933:         ierr    = VecRestoreArray(x,&xarray);
934:         ierr    = VecRestoreArray(y,&yarray);
935:       }
936:     } else if (xnext->type == VECPACK_DA) {
937:       ierr  = VecGetArray(x,&xarray);
938:       ierr  = VecGetArray(y,&yarray);
939:       ierr  = DAGetGlobalVector(xnext->da,&xglobal);
940:       ierr  = DAGetGlobalVector(ynext->da,&yglobal);
941:       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);
942:       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);
943:       if (add) {
944:         ierr  = MatMultAdd(anext->A,xglobal,yglobal,yglobal);
945:       } else {
946:         ierr  = MatMult(anext->A,xglobal,yglobal);
947:       }
948:       ierr  = VecRestoreArray(x,&xarray);
949:       ierr  = VecRestoreArray(y,&yarray);
950:       ierr  = VecResetArray(xglobal);
951:       ierr  = VecResetArray(yglobal);
952:       ierr  = DARestoreGlobalVector(xnext->da,&xglobal);
953:       ierr  = DARestoreGlobalVector(ynext->da,&yglobal);
954:       anext = anext->next;
955:     } else {
956:       SETERRQ(1,"Cannot handle that object type yet");
957:     }
958:     xnext = xnext->next;
959:     ynext = ynext->next;
960:   }
961:   return(0);
962: }

964: int MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
965: {
968:   if (z != y) SETERRQ(1,"Handles y == z only");
969:   MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);
970:   return(0);
971: }

973: int MatMult_Shell_Pack(Mat A,Vec x,Vec y)
974: {
977:   MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);
978:   return(0);
979: }

981: int MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
982: {
983:   struct MatPack     *mpack;
984:   struct VecPackLink *xnext,*ynext;
985:   struct MatPackLink *anext;
986:   Scalar             *xarray,*yarray;
987:   int                ierr;
988:   Vec                xglobal,yglobal;

991:   ierr  = MatShellGetContext(A,(void**)&mpack);
992:   xnext = mpack->left->next;
993:   ynext = mpack->right->next;
994:   anext = mpack->next;

996:   while (xnext) {
997:     if (xnext->type == VECPACK_ARRAY) {
998:       if (!mpack->right->rank) {
999:         ierr    = VecGetArray(x,&xarray);
1000:         ierr    = VecGetArray(y,&yarray);
1001:         ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(Scalar));
1002:         ierr    = VecRestoreArray(x,&xarray);
1003:         ierr    = VecRestoreArray(y,&yarray);
1004:       }
1005:     } else if (xnext->type == VECPACK_DA) {
1006:       ierr  = VecGetArray(x,&xarray);
1007:       ierr  = VecGetArray(y,&yarray);
1008:       ierr  = DAGetGlobalVector(xnext->da,&xglobal);
1009:       ierr  = DAGetGlobalVector(ynext->da,&yglobal);
1010:       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);
1011:       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);
1012:       ierr  = MatMultTranspose(anext->A,xglobal,yglobal);
1013:       ierr  = VecRestoreArray(x,&xarray);
1014:       ierr  = VecRestoreArray(y,&yarray);
1015:       ierr  = VecResetArray(xglobal);
1016:       ierr  = VecResetArray(yglobal);
1017:       ierr  = DARestoreGlobalVector(xnext->da,&xglobal);
1018:       ierr  = DARestoreGlobalVector(ynext->da,&yglobal);
1019:       anext = anext->next;
1020:     } else {
1021:       SETERRQ(1,"Cannot handle that object type yet");
1022:     }
1023:     xnext = xnext->next;
1024:     ynext = ynext->next;
1025:   }
1026:   return(0);
1027: }

1029: int MatDestroy_Shell_Pack(Mat A)
1030: {
1031:   struct MatPack     *mpack;
1032:   struct MatPackLink *anext,*oldanext;
1033:   int                ierr;

1036:   ierr  = MatShellGetContext(A,(void**)&mpack);
1037:   anext = mpack->next;

1039:   while (anext) {
1040:     ierr     = MatDestroy(anext->A);
1041:     oldanext = anext;
1042:     anext    = anext->next;
1043:     ierr     = PetscFree(oldanext);
1044:   }
1045:   PetscFree(mpack);
1046:   return(0);
1047: }

1049: /*@C
1050:     VecPackGetInterpolation - GetInterpolations a VecPack by refining all of its DAs

1052:     Collective on VecPack

1054:     Input Parameters:
1055: +    coarse - coarse grid packer
1056: -    fine - fine grid packer

1058:     Output Parameter:
1059: +    A - interpolation matrix
1060: -    v - scaling vector
1061:  
1062:     Level: advanced

1064: .seealso VecPackDestroy(), VecPackAddArray(), VecPackAddDA(), VecPackCreateGlobalVector(),
1065:          VecPackGather(), VecPackCreate(), VecPackGetGlobalIndices(), VecPackGetAccess()

1067: @*/
1068: int VecPackGetInterpolation(VecPack coarse,VecPack fine,Mat *A,Vec *v)
1069: {
1070:   int                ierr,m,n,M,N;
1071:   struct VecPackLink *nextc  = coarse->next;
1072:   struct VecPackLink *nextf = fine->next;
1073:   struct MatPackLink *nextmat,*pnextmat = 0;
1074:   struct MatPack     *mpack;
1075:   Vec                gcoarse,gfine;

1078:   /* use global vectors only for determining matrix layout */
1079:   VecPackCreateGlobalVector(coarse,&gcoarse);
1080:   VecPackCreateGlobalVector(fine,&gfine);
1081:   VecGetLocalSize(gcoarse,&n);
1082:   VecGetLocalSize(gfine,&m);
1083:   VecGetSize(gcoarse,&N);
1084:   VecGetSize(gfine,&M);
1085:   VecDestroy(gcoarse);
1086:   VecDestroy(gfine);

1088:   ierr         = PetscNew(struct MatPack,&mpack);
1089:   mpack->right = coarse;
1090:   mpack->left  = fine;
1091:   ierr  = MatCreateShell(fine->comm,m,n,M,N,mpack,A);
1092:   ierr  = MatShellSetOperation(*A,MATOP_MULT,(void(*)())MatMult_Shell_Pack);
1093:   ierr  = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)())MatMultTranspose_Shell_Pack);
1094:   ierr  = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)())MatMultAdd_Shell_Pack);
1095:   ierr  = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)())MatDestroy_Shell_Pack);

1097:   /* loop over packed objects, handling one at at time */
1098:   while (nextc) {
1099:     if (nextc->type != nextf->type) SETERRQ(1,"Two VecPack have different layout");

1101:     if (nextc->type == VECPACK_ARRAY) {
1102:       ;
1103:     } else if (nextc->type == VECPACK_DA) {
1104:       ierr          = PetscNew(struct MatPackLink,&nextmat);
1105:       nextmat->next = 0;
1106:       if (pnextmat) {
1107:         pnextmat->next = nextmat;
1108:         pnextmat       = nextmat;
1109:       } else {
1110:         pnextmat    = nextmat;
1111:         mpack->next = nextmat;
1112:       }
1113:       DAGetInterpolation(nextc->da,nextf->da,&nextmat->A,PETSC_NULL);
1114:     } else {
1115:       SETERRQ(1,"Cannot handle that object type yet");
1116:     }
1117:     nextc = nextc->next;
1118:     nextf = nextf->next;
1119:   }
1120:   return(0);
1121: }