Actual source code: shell.c

  2: /*
  3:    This provides a simple shell for Fortran (and C programmers) to 
  4:   create a very simple matrix class for use with KSP without coding 
  5:   much of anything.
  6: */

  8: #include <private/matimpl.h>        /*I "petscmat.h" I*/
  9: #include <private/vecimpl.h>  

 11: typedef struct {
 12:   PetscErrorCode (*destroy)(Mat);
 13:   PetscErrorCode (*mult)(Mat,Vec,Vec);
 14:   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
 15:   PetscErrorCode (*getdiagonal)(Mat,Vec);
 16:   PetscScalar    vscale,vshift;
 17:   Vec            dshift;
 18:   Vec            left,right;
 19:   Vec            dshift_owned,left_owned,right_owned;
 20:   Vec            left_work,right_work;
 21:   PetscBool      usingscaled;
 22:   void           *ctx;
 23: } Mat_Shell;
 24: /*
 25:  The most general expression for the matrix is

 27:  S = L (a A + B) R

 29:  where
 30:  A is the matrix defined by the user's function
 31:  a is a scalar multiple
 32:  L is left scaling
 33:  R is right scaling
 34:  B is a diagonal shift defined by
 35:    diag(dshift) if the vector dshift is non-NULL
 36:    vscale*identity otherwise

 38:  The following identities apply:

 40:  Scale by c:
 41:   c [L (a A + B) R] = L [(a c) A + c B] R

 43:  Shift by c:
 44:   [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R

 46:  Diagonal scaling is achieved by simply multiplying with existing L and R vectors

 48:  In the data structure:

 50:  vscale=1.0  means no special scaling will be applied
 51:  dshift=NULL means a constant diagonal shift (fall back to vshift)
 52:  vshift=0.0  means no constant diagonal shift, note that vshift is only used if dshift is NULL
 53: */

 55: static PetscErrorCode MatMult_Shell(Mat,Vec,Vec);
 56: static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec);
 57: static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec);

 61: static PetscErrorCode MatShellUseScaledMethods(Mat Y)
 62: {
 63:   Mat_Shell *shell = (Mat_Shell*)Y->data;

 66:   if (shell->usingscaled) return(0);
 67:   shell->mult  = Y->ops->mult;
 68:   Y->ops->mult = MatMult_Shell;
 69:   if (Y->ops->multtranspose) {
 70:     shell->multtranspose  = Y->ops->multtranspose;
 71:     Y->ops->multtranspose = MatMultTranspose_Shell;
 72:   }
 73:   if (Y->ops->getdiagonal) {
 74:     shell->getdiagonal  = Y->ops->getdiagonal;
 75:     Y->ops->getdiagonal = MatGetDiagonal_Shell;
 76:   }
 77:   shell->usingscaled = PETSC_TRUE;
 78:   return(0);
 79: }

 83: static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
 84: {
 85:   Mat_Shell *shell = (Mat_Shell*)A->data;

 89:   *xx = PETSC_NULL;
 90:   if (!shell->left) {
 91:     *xx = x;
 92:   } else {
 93:     if (!shell->left_work) {VecDuplicate(shell->left,&shell->left_work);}
 94:     VecPointwiseMult(shell->left_work,x,shell->left);
 95:     *xx = shell->left_work;
 96:   }
 97:   return(0);
 98: }

102: static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
103: {
104:   Mat_Shell *shell = (Mat_Shell*)A->data;

108:   *xx = PETSC_NULL;
109:   if (!shell->right) {
110:     *xx = x;
111:   } else {
112:     if (!shell->right_work) {VecDuplicate(shell->right,&shell->right_work);}
113:     VecPointwiseMult(shell->right_work,x,shell->right);
114:     *xx = shell->right_work;
115:   }
116:   return(0);
117: }

121: static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
122: {
123:   Mat_Shell *shell = (Mat_Shell*)A->data;

127:   if (shell->left) {VecPointwiseMult(x,x,shell->left);}
128:   return(0);
129: }

133: static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
134: {
135:   Mat_Shell *shell = (Mat_Shell*)A->data;

139:   if (shell->right) {VecPointwiseMult(x,x,shell->right);}
140:   return(0);
141: }

145: static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
146: {
147:   Mat_Shell *shell = (Mat_Shell*)A->data;

151:   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
152:     PetscInt          i,m;
153:     const PetscScalar *x,*d;
154:     PetscScalar       *y;
155:     VecGetLocalSize(X,&m);
156:     VecGetArrayRead(shell->dshift,&d);
157:     VecGetArrayRead(X,&x);
158:     VecGetArray(Y,&y);
159:     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
160:     VecRestoreArrayRead(shell->dshift,&d);
161:     VecRestoreArrayRead(X,&x);
162:     VecRestoreArray(Y,&y);
163:   } else if (PetscAbsScalar(shell->vshift) != 0) {
164:     VecAXPBY(Y,shell->vshift,shell->vscale,X);
165:   } else if (shell->vscale != 1.0) {
166:     VecScale(Y,shell->vscale);
167:   }
168:   return(0);
169: }

173: /*@C
174:     MatShellGetContext - Returns the user-provided context associated with a shell matrix.

176:     Not Collective

178:     Input Parameter:
179: .   mat - the matrix, should have been created with MatCreateShell()

181:     Output Parameter:
182: .   ctx - the user provided context

184:     Level: advanced

186:     Notes:
187:     This routine is intended for use within various shell matrix routines,
188:     as set with MatShellSetOperation().
189:     
190: .keywords: matrix, shell, get, context

192: .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
193: @*/
194: PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
195: {
197:   PetscBool      flg;

202:   PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);
203:   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
204:   else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
205:   return(0);
206: }

210: PetscErrorCode MatDestroy_Shell(Mat mat)
211: {
213:   Mat_Shell      *shell = (Mat_Shell*)mat->data;

216:   if (shell->destroy) {
217:     (*shell->destroy)(mat);
218:   }
219:   VecDestroy(&shell->left_owned);
220:   VecDestroy(&shell->right_owned);
221:   VecDestroy(&shell->dshift_owned);
222:   VecDestroy(&shell->left_work);
223:   VecDestroy(&shell->right_work);
224:   PetscFree(mat->data);
225:   return(0);
226: }

230: PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
231: {
232:   Mat_Shell      *shell = (Mat_Shell*)A->data;
234:   Vec            xx;

237:   MatShellPreScaleRight(A,x,&xx);
238:   (*shell->mult)(A,xx,y);
239:   MatShellShiftAndScale(A,xx,y);
240:   MatShellPostScaleLeft(A,y);
241:   return(0);
242: }

246: PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
247: {
248:   Mat_Shell      *shell = (Mat_Shell*)A->data;
250:   Vec            xx;

253:   MatShellPreScaleLeft(A,x,&xx);
254:   (*shell->multtranspose)(A,xx,y);
255:   MatShellShiftAndScale(A,xx,y);
256:   MatShellPostScaleRight(A,y);
257:   return(0);
258: }

262: PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
263: {
264:   Mat_Shell      *shell = (Mat_Shell*)A->data;

268:   (*shell->getdiagonal)(A,v);
269:   VecScale(v,shell->vscale);
270:   if (shell->dshift) {
271:     VecPointwiseMult(v,v,shell->dshift);
272:   } else {
273:     VecShift(v,shell->vshift);
274:   }
275:   if (shell->left)  {VecPointwiseMult(v,v,shell->left);}
276:   if (shell->right) {VecPointwiseMult(v,v,shell->right);}
277:   return(0);
278: }

282: PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
283: {
284:   Mat_Shell *shell = (Mat_Shell*)Y->data;

288:   if (shell->left || shell->right || shell->dshift) {
289:     if (!shell->dshift) {
290:       if (!shell->dshift_owned) {VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);}
291:       shell->dshift = shell->dshift_owned;
292:       VecSet(shell->dshift,shell->vshift+a);
293:     } else {VecScale(shell->dshift,a);}
294:     if (shell->left)  {VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);}
295:     if (shell->right) {VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);}
296:   } else shell->vshift += a;
297:   MatShellUseScaledMethods(Y);
298:   return(0);
299: }

303: PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
304: {
305:   Mat_Shell      *shell = (Mat_Shell*)Y->data;

309:   shell->vscale *= a;
310:   if (shell->dshift) {VecScale(shell->dshift,a);}
311:   else shell->vshift *= a;
312:   MatShellUseScaledMethods(Y);
313:   return(0);
314: }

318: static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
319: {
320:   Mat_Shell *shell = (Mat_Shell*)Y->data;

324:   if (left) {
325:     if (!shell->left_owned) {VecDuplicate(left,&shell->left_owned);}
326:     if (shell->left) {VecPointwiseMult(shell->left,shell->left,left);}
327:     else {
328:       shell->left = shell->left_owned;
329:       VecCopy(left,shell->left);
330:     }
331:   }
332:   if (right) {
333:     if (!shell->right_owned) {VecDuplicate(right,&shell->right_owned);}
334:     if (shell->right) {VecPointwiseMult(shell->right,shell->right,right);}
335:     else {
336:       shell->right = shell->right_owned;
337:       VecCopy(right,shell->right);
338:     }
339:   }
340:   MatShellUseScaledMethods(Y);
341:   return(0);
342: }

346: PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
347: {
348:   Mat_Shell *shell = (Mat_Shell*)Y->data;

351:   if (t == MAT_FINAL_ASSEMBLY) {
352:     shell->vshift      = 0.0;
353:     shell->vscale      = 1.0;
354:     shell->dshift      = PETSC_NULL;
355:     shell->left        = PETSC_NULL;
356:     shell->right       = PETSC_NULL;
357:     if (shell->mult) {
358:       Y->ops->mult = shell->mult;
359:       shell->mult  = PETSC_NULL;
360:     }
361:     if (shell->multtranspose) {
362:       Y->ops->multtranspose = shell->multtranspose;
363:       shell->multtranspose  = PETSC_NULL;
364:     }
365:     if (shell->getdiagonal) {
366:       Y->ops->getdiagonal = shell->getdiagonal;
367:       shell->getdiagonal  = PETSC_NULL;
368:     }
369:     shell->usingscaled = PETSC_FALSE;
370:   }
371:   return(0);
372: }


378: PetscErrorCode MatSetBlockSize_Shell(Mat A,PetscInt bs)
379: {

383:   PetscLayoutSetBlockSize(A->rmap,bs);
384:   PetscLayoutSetBlockSize(A->cmap,bs);
385:   return(0);
386: }

388: static struct _MatOps MatOps_Values = {0,
389:        0,
390:        0,
391:        0,
392: /* 4*/ 0,
393:        0,
394:        0,
395:        0,
396:        0,
397:        0,
398: /*10*/ 0,
399:        0,
400:        0,
401:        0,
402:        0,
403: /*15*/ 0,
404:        0,
405:        0,
406:        MatDiagonalScale_Shell,
407:        0,
408: /*20*/ 0,
409:        MatAssemblyEnd_Shell,
410:        0,
411:        0,
412: /*24*/ 0,
413:        0,
414:        0,
415:        0,
416:        0,
417: /*29*/ 0,
418:        0,
419:        0,
420:        0,
421:        0,
422: /*34*/ 0,
423:        0,
424:        0,
425:        0,
426:        0,
427: /*39*/ 0,
428:        0,
429:        0,
430:        0,
431:        0,
432: /*44*/ 0,
433:        MatScale_Shell,
434:        MatShift_Shell,
435:        0,
436:        0,
437: /*49*/ MatSetBlockSize_Shell,
438:        0,
439:        0,
440:        0,
441:        0,
442: /*54*/ 0,
443:        0,
444:        0,
445:        0,
446:        0,
447: /*59*/ 0,
448:        MatDestroy_Shell,
449:        0,
450:        0,
451:        0,
452: /*64*/ 0,
453:        0,
454:        0,
455:        0,
456:        0,
457: /*69*/ 0,
458:        0,
459:        MatConvert_Shell,
460:        0,
461:        0,
462: /*74*/ 0,
463:        0,
464:        0,
465:        0,
466:        0,
467: /*79*/ 0,
468:        0,
469:        0,
470:        0,
471:        0,
472: /*84*/ 0,
473:        0,
474:        0,
475:        0,
476:        0,
477: /*89*/ 0,
478:        0,
479:        0,
480:        0,
481:        0,
482: /*94*/ 0,
483:        0,
484:        0,
485:        0};

487: /*MC
488:    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.

490:   Level: advanced

492: .seealso: MatCreateShell
493: M*/

498: PetscErrorCode  MatCreate_Shell(Mat A)
499: {
500:   Mat_Shell      *b;

504:   PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));

506:   PetscNewLog(A,Mat_Shell,&b);
507:   A->data = (void*)b;

509:   PetscLayoutSetBlockSize(A->rmap,1);
510:   PetscLayoutSetBlockSize(A->cmap,1);
511:   PetscLayoutSetUp(A->rmap);
512:   PetscLayoutSetUp(A->cmap);

514:   b->ctx           = 0;
515:   b->vshift        = 0.0;
516:   b->vscale        = 1.0;
517:   b->mult          = 0;
518:   b->multtranspose = 0;
519:   b->getdiagonal   = 0;
520:   A->assembled     = PETSC_TRUE;
521:   A->preallocated  = PETSC_FALSE;
522:   PetscObjectChangeTypeName((PetscObject)A,MATSHELL);
523:   return(0);
524: }

529: /*@C
530:    MatCreateShell - Creates a new matrix class for use with a user-defined
531:    private data storage format. 

533:   Collective on MPI_Comm

535:    Input Parameters:
536: +  comm - MPI communicator
537: .  m - number of local rows (must be given)
538: .  n - number of local columns (must be given)
539: .  M - number of global rows (may be PETSC_DETERMINE)
540: .  N - number of global columns (may be PETSC_DETERMINE)
541: -  ctx - pointer to data needed by the shell matrix routines

543:    Output Parameter:
544: .  A - the matrix

546:    Level: advanced

548:   Usage:
550: $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
551: $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
552: $    [ Use matrix for operations that have been set ]
553: $    MatDestroy(mat);

555:    Notes:
556:    The shell matrix type is intended to provide a simple class to use
557:    with KSP (such as, for use with matrix-free methods). You should not
558:    use the shell type if you plan to define a complete matrix class.

560:    Fortran Notes: The context can only be an integer or a PetscObject
561:       unfortunately it cannot be a Fortran array or derived type.

563:    PETSc requires that matrices and vectors being used for certain
564:    operations are partitioned accordingly.  For example, when
565:    creating a shell matrix, A, that supports parallel matrix-vector
566:    products using MatMult(A,x,y) the user should set the number
567:    of local matrix rows to be the number of local elements of the
568:    corresponding result vector, y. Note that this is information is
569:    required for use of the matrix interface routines, even though
570:    the shell matrix may not actually be physically partitioned.
571:    For example,

573: $
574: $     Vec x, y
576: $     Mat A
577: $
578: $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
579: $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
580: $     VecGetLocalSize(y,&m);
581: $     VecGetLocalSize(x,&n);
582: $     MatCreateShell(comm,m,n,M,N,ctx,&A);
583: $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
584: $     MatMult(A,x,y);
585: $     MatDestroy(A);
586: $     VecDestroy(y); VecDestroy(x);
587: $

589: .keywords: matrix, shell, create

591: .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
592: @*/
593: PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
594: {

598:   MatCreate(comm,A);
599:   MatSetSizes(*A,m,n,M,N);
600: 
601:   MatSetType(*A,MATSHELL);
602:   MatShellSetContext(*A,ctx);
603:   return(0);
604: }

608: /*@
609:     MatShellSetContext - sets the context for a shell matrix

611:    Logically Collective on Mat

613:     Input Parameters:
614: +   mat - the shell matrix
615: -   ctx - the context

617:    Level: advanced

619:    Fortran Notes: The context can only be an integer or a PetscObject
620:       unfortunately it cannot be a Fortran array or derived type.

622: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
623: @*/
624: PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
625: {
626:   Mat_Shell      *shell = (Mat_Shell*)mat->data;
628:   PetscBool      flg;

632:   PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);
633:   if (flg) {
634:     shell->ctx = ctx;
635:   } else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
636:   return(0);
637: }

641: /*@C
642:     MatShellSetOperation - Allows user to set a matrix operation for
643:                            a shell matrix.

645:    Logically Collective on Mat

647:     Input Parameters:
648: +   mat - the shell matrix
649: .   op - the name of the operation
650: -   f - the function that provides the operation.

652:    Level: advanced

654:     Usage:
656: $      MatCreateShell(comm,m,n,M,N,ctx,&A);
657: $      MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);

659:     Notes:
660:     See the file include/petscmat.h for a complete list of matrix
661:     operations, which all have the form MATOP_<OPERATION>, where
662:     <OPERATION> is the name (in all capital letters) of the
663:     user interface routine (e.g., MatMult() -> MATOP_MULT).

665:     All user-provided functions should have the same calling
666:     sequence as the usual matrix interface routines, since they
667:     are intended to be accessed via the usual matrix interface
668:     routines, e.g., 
669: $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)

671:     In particular each function MUST return an error code of 0 on success and 
672:     nonzero on failure.

674:     Within each user-defined routine, the user should call
675:     MatShellGetContext() to obtain the user-defined context that was
676:     set by MatCreateShell().

678:     Fortran Notes: For MatGetVecs() the user code should check if the input left or right matrix is -1 and in that case not
679:        generate a matrix. See src/mat/examples/tests/ex120f.F

681: .keywords: matrix, shell, set, operation

683: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
684: @*/
685: PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
686: {
688:   PetscBool      flg;

692:   if (op == MATOP_DESTROY) {
693:     PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);
694:     if (flg) {
695:        Mat_Shell *shell = (Mat_Shell*)mat->data;
696:        shell->destroy                 = (PetscErrorCode (*)(Mat)) f;
697:     } else mat->ops->destroy          = (PetscErrorCode (*)(Mat)) f;
698:   }
699:   else if (op == MATOP_VIEW) mat->ops->view  = (PetscErrorCode (*)(Mat,PetscViewer)) f;
700:   else                       (((void(**)(void))mat->ops)[op]) = f;

702:   return(0);
703: }

707: /*@C
708:     MatShellGetOperation - Gets a matrix function for a shell matrix.

710:     Not Collective

712:     Input Parameters:
713: +   mat - the shell matrix
714: -   op - the name of the operation

716:     Output Parameter:
717: .   f - the function that provides the operation.

719:     Level: advanced

721:     Notes:
722:     See the file include/petscmat.h for a complete list of matrix
723:     operations, which all have the form MATOP_<OPERATION>, where
724:     <OPERATION> is the name (in all capital letters) of the
725:     user interface routine (e.g., MatMult() -> MATOP_MULT).

727:     All user-provided functions have the same calling
728:     sequence as the usual matrix interface routines, since they
729:     are intended to be accessed via the usual matrix interface
730:     routines, e.g., 
731: $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)

733:     Within each user-defined routine, the user should call
734:     MatShellGetContext() to obtain the user-defined context that was
735:     set by MatCreateShell().

737: .keywords: matrix, shell, set, operation

739: .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
740: @*/
741: PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
742: {
744:   PetscBool      flg;

748:   if (op == MATOP_DESTROY) {
749:     PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);
750:     if (flg) {
751:       Mat_Shell *shell = (Mat_Shell*)mat->data;
752:       *f = (void(*)(void))shell->destroy;
753:     } else {
754:       *f = (void(*)(void))mat->ops->destroy;
755:     }
756:   } else if (op == MATOP_VIEW) {
757:     *f = (void(*)(void))mat->ops->view;
758:   } else {
759:     *f = (((void(**)(void))mat->ops)[op]);
760:   }

762:   return(0);
763: }