xref: /petsc/src/mat/impls/shell/shell.c (revision cb8c8a77bf3139dca4664c3836b73ceb6d6ddabe)
1be1d678aSKris Buschelman 
2e51e0e81SBarry Smith /*
320563c6bSBarry Smith    This provides a simple shell for Fortran (and C programmers) to
420563c6bSBarry Smith   create a very simple matrix class for use with KSP without coding
5ed3cc1f0SBarry Smith   much of anything.
6e51e0e81SBarry Smith */
7e51e0e81SBarry Smith 
8af0996ceSBarry Smith #include <petsc/private/matimpl.h>        /*I "petscmat.h" I*/
9e51e0e81SBarry Smith 
1028f357e3SAlex Fikl struct _MatShellOps {
1128f357e3SAlex Fikl   /*   0 */
126849ba73SBarry Smith   PetscErrorCode (*mult)(Mat,Vec,Vec);
1328f357e3SAlex Fikl   /*   5 */
1418be62a5SSatish Balay   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
1528f357e3SAlex Fikl   /*  10 */
1628f357e3SAlex Fikl   /*  15 */
1718be62a5SSatish Balay   PetscErrorCode (*getdiagonal)(Mat,Vec);
1828f357e3SAlex Fikl   /*  20 */
1928f357e3SAlex Fikl   /*  24 */
2028f357e3SAlex Fikl   /*  29 */
2128f357e3SAlex Fikl   /*  34 */
2228f357e3SAlex Fikl   /*  39 */
2328f357e3SAlex Fikl   PetscErrorCode (*copy)(Mat,Mat,MatStructure);
2428f357e3SAlex Fikl   /*  44 */
256464bf51SAlex Fikl   PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode);
2628f357e3SAlex Fikl   /*  49 */
2728f357e3SAlex Fikl   /*  54 */
2828f357e3SAlex Fikl   /*  59 */
2928f357e3SAlex Fikl   PetscErrorCode (*destroy)(Mat);
3028f357e3SAlex Fikl   /*  64 */
3128f357e3SAlex Fikl   /*  69 */
3228f357e3SAlex Fikl   /*  74 */
3328f357e3SAlex Fikl   /*  79 */
3428f357e3SAlex Fikl   /*  84 */
3528f357e3SAlex Fikl   /*  89 */
3628f357e3SAlex Fikl   /*  94 */
3728f357e3SAlex Fikl   /*  99 */
3828f357e3SAlex Fikl   /* 104 */
3928f357e3SAlex Fikl   /* 109 */
4028f357e3SAlex Fikl   /* 114 */
4128f357e3SAlex Fikl   /* 119 */
4228f357e3SAlex Fikl   /* 124 */
4328f357e3SAlex Fikl   /* 129 */
4428f357e3SAlex Fikl   /* 134 */
4528f357e3SAlex Fikl   /* 139 */
4628f357e3SAlex Fikl   /* 144 */
4728f357e3SAlex Fikl };
4828f357e3SAlex Fikl 
4928f357e3SAlex Fikl typedef struct {
5028f357e3SAlex Fikl   struct _MatShellOps ops[1];
512205254eSKarl Rupp 
52ef66eb69SBarry Smith   PetscScalar vscale,vshift;
538fe8eb27SJed Brown   Vec         dshift;
548fe8eb27SJed Brown   Vec         left,right;
558fe8eb27SJed Brown   Vec         dshift_owned,left_owned,right_owned;
568fe8eb27SJed Brown   Vec         left_work,right_work;
575edf6516SJed Brown   Vec         left_add_work,right_add_work;
588fe8eb27SJed Brown   PetscBool   usingscaled;
5920563c6bSBarry Smith   void        *ctx;
6088cf3e7dSBarry Smith } Mat_Shell;
618fe8eb27SJed Brown /*
628fe8eb27SJed Brown  The most general expression for the matrix is
638fe8eb27SJed Brown 
648fe8eb27SJed Brown  S = L (a A + B) R
658fe8eb27SJed Brown 
668fe8eb27SJed Brown  where
678fe8eb27SJed Brown  A is the matrix defined by the user's function
688fe8eb27SJed Brown  a is a scalar multiple
698fe8eb27SJed Brown  L is left scaling
708fe8eb27SJed Brown  R is right scaling
718fe8eb27SJed Brown  B is a diagonal shift defined by
728fe8eb27SJed Brown    diag(dshift) if the vector dshift is non-NULL
738fe8eb27SJed Brown    vscale*identity otherwise
748fe8eb27SJed Brown 
758fe8eb27SJed Brown  The following identities apply:
768fe8eb27SJed Brown 
778fe8eb27SJed Brown  Scale by c:
788fe8eb27SJed Brown   c [L (a A + B) R] = L [(a c) A + c B] R
798fe8eb27SJed Brown 
808fe8eb27SJed Brown  Shift by c:
818fe8eb27SJed Brown   [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R
828fe8eb27SJed Brown 
838fe8eb27SJed Brown  Diagonal scaling is achieved by simply multiplying with existing L and R vectors
848fe8eb27SJed Brown 
858fe8eb27SJed Brown  In the data structure:
868fe8eb27SJed Brown 
878fe8eb27SJed Brown  vscale=1.0  means no special scaling will be applied
888fe8eb27SJed Brown  dshift=NULL means a constant diagonal shift (fall back to vshift)
898fe8eb27SJed Brown  vshift=0.0  means no constant diagonal shift, note that vshift is only used if dshift is NULL
908fe8eb27SJed Brown */
918fe8eb27SJed Brown 
928fe8eb27SJed Brown static PetscErrorCode MatMult_Shell(Mat,Vec,Vec);
938fe8eb27SJed Brown static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec);
948fe8eb27SJed Brown static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec);
9528f357e3SAlex Fikl static PetscErrorCode MatCopy_Shell(Mat,Mat,MatStructure);
968fe8eb27SJed Brown 
978fe8eb27SJed Brown static PetscErrorCode MatShellUseScaledMethods(Mat Y)
988fe8eb27SJed Brown {
998fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)Y->data;
1008fe8eb27SJed Brown 
1018fe8eb27SJed Brown   PetscFunctionBegin;
1028fe8eb27SJed Brown   if (shell->usingscaled) PetscFunctionReturn(0);
10328f357e3SAlex Fikl   shell->ops->mult  = Y->ops->mult;
1048fe8eb27SJed Brown   Y->ops->mult = MatMult_Shell;
1058fe8eb27SJed Brown   if (Y->ops->multtranspose) {
10628f357e3SAlex Fikl     shell->ops->multtranspose  = Y->ops->multtranspose;
1078fe8eb27SJed Brown     Y->ops->multtranspose = MatMultTranspose_Shell;
1088fe8eb27SJed Brown   }
1098fe8eb27SJed Brown   if (Y->ops->getdiagonal) {
11028f357e3SAlex Fikl     shell->ops->getdiagonal  = Y->ops->getdiagonal;
1118fe8eb27SJed Brown     Y->ops->getdiagonal = MatGetDiagonal_Shell;
1128fe8eb27SJed Brown   }
11328f357e3SAlex Fikl   if (Y->ops->copy) {
11428f357e3SAlex Fikl     shell->ops->copy  = Y->ops->copy;
11528f357e3SAlex Fikl     Y->ops->copy = MatCopy_Shell;
11628f357e3SAlex Fikl   }
1178fe8eb27SJed Brown   shell->usingscaled = PETSC_TRUE;
1188fe8eb27SJed Brown   PetscFunctionReturn(0);
1198fe8eb27SJed Brown }
1208fe8eb27SJed Brown 
1218fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
1228fe8eb27SJed Brown {
1238fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1248fe8eb27SJed Brown   PetscErrorCode ierr;
1258fe8eb27SJed Brown 
1268fe8eb27SJed Brown   PetscFunctionBegin;
1270298fd71SBarry Smith   *xx = NULL;
1288fe8eb27SJed Brown   if (!shell->left) {
1298fe8eb27SJed Brown     *xx = x;
1308fe8eb27SJed Brown   } else {
1318fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
1328fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
1338fe8eb27SJed Brown     *xx  = shell->left_work;
1348fe8eb27SJed Brown   }
1358fe8eb27SJed Brown   PetscFunctionReturn(0);
1368fe8eb27SJed Brown }
1378fe8eb27SJed Brown 
1388fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
1398fe8eb27SJed Brown {
1408fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1418fe8eb27SJed Brown   PetscErrorCode ierr;
1428fe8eb27SJed Brown 
1438fe8eb27SJed Brown   PetscFunctionBegin;
1440298fd71SBarry Smith   *xx = NULL;
1458fe8eb27SJed Brown   if (!shell->right) {
1468fe8eb27SJed Brown     *xx = x;
1478fe8eb27SJed Brown   } else {
1488fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
1498fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
1508fe8eb27SJed Brown     *xx  = shell->right_work;
1518fe8eb27SJed Brown   }
1528fe8eb27SJed Brown   PetscFunctionReturn(0);
1538fe8eb27SJed Brown }
1548fe8eb27SJed Brown 
1558fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
1568fe8eb27SJed Brown {
1578fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1588fe8eb27SJed Brown   PetscErrorCode ierr;
1598fe8eb27SJed Brown 
1608fe8eb27SJed Brown   PetscFunctionBegin;
1618fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
1628fe8eb27SJed Brown   PetscFunctionReturn(0);
1638fe8eb27SJed Brown }
1648fe8eb27SJed Brown 
1658fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
1668fe8eb27SJed Brown {
1678fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1688fe8eb27SJed Brown   PetscErrorCode ierr;
1698fe8eb27SJed Brown 
1708fe8eb27SJed Brown   PetscFunctionBegin;
1718fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
1728fe8eb27SJed Brown   PetscFunctionReturn(0);
1738fe8eb27SJed Brown }
1748fe8eb27SJed Brown 
1758fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
1768fe8eb27SJed Brown {
1778fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1788fe8eb27SJed Brown   PetscErrorCode ierr;
1798fe8eb27SJed Brown 
1808fe8eb27SJed Brown   PetscFunctionBegin;
1818fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
1828fe8eb27SJed Brown     PetscInt          i,m;
1838fe8eb27SJed Brown     const PetscScalar *x,*d;
1848fe8eb27SJed Brown     PetscScalar       *y;
1858fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
1868fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1878fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1888fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1898fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
1908fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1918fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1928fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
193880538faSHong Zhang   } else if (PetscAbsScalar(shell->vshift) != 0) {
1948fe8eb27SJed Brown     ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr);
195027c5fb5SJed Brown   } else if (shell->vscale != 1.0) {
196027c5fb5SJed Brown     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
1978fe8eb27SJed Brown   }
1988fe8eb27SJed Brown   PetscFunctionReturn(0);
1998fe8eb27SJed Brown }
200e51e0e81SBarry Smith 
2019d225801SJed Brown /*@
202a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
203b4fd4287SBarry Smith 
204c7fcc2eaSBarry Smith     Not Collective
205c7fcc2eaSBarry Smith 
206b4fd4287SBarry Smith     Input Parameter:
207b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
208b4fd4287SBarry Smith 
209b4fd4287SBarry Smith     Output Parameter:
210b4fd4287SBarry Smith .   ctx - the user provided context
211b4fd4287SBarry Smith 
21215091d37SBarry Smith     Level: advanced
21315091d37SBarry Smith 
214daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
215daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
216a62d957aSLois Curfman McInnes 
217a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context
218a62d957aSLois Curfman McInnes 
219ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
2200bc0a719SBarry Smith @*/
2218fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
222b4fd4287SBarry Smith {
223dfbe8321SBarry Smith   PetscErrorCode ierr;
224ace3abfcSBarry Smith   PetscBool      flg;
225273d9f13SBarry Smith 
2263a40ed3dSBarry Smith   PetscFunctionBegin;
2270700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2284482741eSBarry Smith   PetscValidPointer(ctx,2);
229251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
230940183f3SBarry Smith   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
231ce94432eSBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
2323a40ed3dSBarry Smith   PetscFunctionReturn(0);
233b4fd4287SBarry Smith }
234b4fd4287SBarry Smith 
235dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
236e51e0e81SBarry Smith {
237dfbe8321SBarry Smith   PetscErrorCode ierr;
238bf0cc555SLisandro Dalcin   Mat_Shell      *shell = (Mat_Shell*)mat->data;
239ed3cc1f0SBarry Smith 
2403a40ed3dSBarry Smith   PetscFunctionBegin;
24128f357e3SAlex Fikl   if (shell->ops->destroy) {
24228f357e3SAlex Fikl     ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr);
243bf0cc555SLisandro Dalcin   }
2448fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr);
2458fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr);
2468fe8eb27SJed Brown   ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr);
2478fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
2488fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
2495edf6516SJed Brown   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
2505edf6516SJed Brown   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
251bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
2523a40ed3dSBarry Smith   PetscFunctionReturn(0);
253e51e0e81SBarry Smith }
254e51e0e81SBarry Smith 
2557fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
2567fabda3fSAlex Fikl {
25728f357e3SAlex Fikl   Mat_Shell       *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
2587fabda3fSAlex Fikl   PetscErrorCode  ierr;
259*cb8c8a77SBarry Smith   PetscBool       matflg;
2607fabda3fSAlex Fikl 
2617fabda3fSAlex Fikl   PetscFunctionBegin;
26228f357e3SAlex Fikl   ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr);
263*cb8c8a77SBarry Smith   if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL");
2647fabda3fSAlex Fikl 
2657fabda3fSAlex Fikl   ierr = MatShellUseScaledMethods(B);CHKERRQ(ierr);
266*cb8c8a77SBarry Smith   ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
267*cb8c8a77SBarry Smith   ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
2687fabda3fSAlex Fikl 
269*cb8c8a77SBarry Smith   if (shellA->ops->copy) {
27028f357e3SAlex Fikl     ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr);
271*cb8c8a77SBarry Smith   }
2727fabda3fSAlex Fikl   shellB->vscale = shellA->vscale;
2737fabda3fSAlex Fikl   shellB->vshift = shellA->vshift;
2747fabda3fSAlex Fikl   if (shellA->dshift_owned) {
2757fabda3fSAlex Fikl     if (!shellB->dshift_owned) {
2767fabda3fSAlex Fikl       ierr = VecDuplicate(shellA->dshift_owned,&shellB->dshift_owned);CHKERRQ(ierr);
2777fabda3fSAlex Fikl     }
2787fabda3fSAlex Fikl     ierr = VecCopy(shellA->dshift_owned,shellB->dshift_owned);CHKERRQ(ierr);
2797fabda3fSAlex Fikl     shellB->dshift = shellB->dshift_owned;
2807fabda3fSAlex Fikl   } else {
2817fabda3fSAlex Fikl     shellB->dshift = NULL;
2827fabda3fSAlex Fikl   }
2837fabda3fSAlex Fikl   if (shellA->left_owned) {
2847fabda3fSAlex Fikl     if (!shellB->left_owned) {
2857fabda3fSAlex Fikl       ierr = VecDuplicate(shellA->left_owned,&shellB->left_owned);CHKERRQ(ierr);
2867fabda3fSAlex Fikl     }
2877fabda3fSAlex Fikl     ierr = VecCopy(shellA->left_owned,shellB->left_owned);CHKERRQ(ierr);
2887fabda3fSAlex Fikl     shellB->left = shellB->left_owned;
2897fabda3fSAlex Fikl   } else {
2907fabda3fSAlex Fikl     shellB->left = NULL;
2917fabda3fSAlex Fikl   }
2927fabda3fSAlex Fikl   if (shellA->right_owned) {
2937fabda3fSAlex Fikl     if (!shellB->right_owned) {
2947fabda3fSAlex Fikl       ierr = VecDuplicate(shellA->right_owned,&shellB->right_owned);CHKERRQ(ierr);
2957fabda3fSAlex Fikl     }
2967fabda3fSAlex Fikl     ierr = VecCopy(shellA->right_owned,shellB->right_owned);CHKERRQ(ierr);
2977fabda3fSAlex Fikl     shellB->right = shellB->right_owned;
2987fabda3fSAlex Fikl   } else {
2997fabda3fSAlex Fikl     shellB->right = NULL;
3007fabda3fSAlex Fikl   }
3017fabda3fSAlex Fikl   PetscFunctionReturn(0);
3027fabda3fSAlex Fikl }
3037fabda3fSAlex Fikl 
304*cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M)
305*cb8c8a77SBarry Smith {
306*cb8c8a77SBarry Smith   PetscErrorCode ierr;
307*cb8c8a77SBarry Smith   void           *ctx;
308*cb8c8a77SBarry Smith 
309*cb8c8a77SBarry Smith   PetscFunctionBegin;
310*cb8c8a77SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
311*cb8c8a77SBarry Smith   ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr);
312*cb8c8a77SBarry Smith   ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
313*cb8c8a77SBarry Smith   PetscFunctionReturn(0);
314*cb8c8a77SBarry Smith }
315*cb8c8a77SBarry Smith 
316dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
317ef66eb69SBarry Smith {
318ef66eb69SBarry Smith   Mat_Shell        *shell = (Mat_Shell*)A->data;
319dfbe8321SBarry Smith   PetscErrorCode   ierr;
32025578ef6SJed Brown   Vec              xx;
321e3f487b0SBarry Smith   PetscObjectState instate,outstate;
322ef66eb69SBarry Smith 
323ef66eb69SBarry Smith   PetscFunctionBegin;
3248fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
325e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
32628f357e3SAlex Fikl   ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr);
327e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
328e3f487b0SBarry Smith   if (instate == outstate) {
329e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
330e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
331e3f487b0SBarry Smith   }
3328fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
3338fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
334ef66eb69SBarry Smith   PetscFunctionReturn(0);
335ef66eb69SBarry Smith }
336ef66eb69SBarry Smith 
3375edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3385edf6516SJed Brown {
3395edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3405edf6516SJed Brown   PetscErrorCode ierr;
3415edf6516SJed Brown 
3425edf6516SJed Brown   PetscFunctionBegin;
3435edf6516SJed Brown   if (y == z) {
3445edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
3455edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
346b8430d49SBarry Smith     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
3475edf6516SJed Brown   } else {
3485edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
3495edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3505edf6516SJed Brown   }
3515edf6516SJed Brown   PetscFunctionReturn(0);
3525edf6516SJed Brown }
3535edf6516SJed Brown 
35418be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
35518be62a5SSatish Balay {
35618be62a5SSatish Balay   Mat_Shell        *shell = (Mat_Shell*)A->data;
35718be62a5SSatish Balay   PetscErrorCode   ierr;
35825578ef6SJed Brown   Vec              xx;
359e3f487b0SBarry Smith   PetscObjectState instate,outstate;
36018be62a5SSatish Balay 
36118be62a5SSatish Balay   PetscFunctionBegin;
3628fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
363e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
36428f357e3SAlex Fikl   ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr);
365e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
366e3f487b0SBarry Smith   if (instate == outstate) {
367e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
368e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
369e3f487b0SBarry Smith   }
3708fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
3718fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
37218be62a5SSatish Balay   PetscFunctionReturn(0);
37318be62a5SSatish Balay }
37418be62a5SSatish Balay 
3755edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3765edf6516SJed Brown {
3775edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3785edf6516SJed Brown   PetscErrorCode ierr;
3795edf6516SJed Brown 
3805edf6516SJed Brown   PetscFunctionBegin;
3815edf6516SJed Brown   if (y == z) {
3825edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
3835edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
3845edf6516SJed Brown     ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr);
3855edf6516SJed Brown   } else {
3865edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
3875edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3885edf6516SJed Brown   }
3895edf6516SJed Brown   PetscFunctionReturn(0);
3905edf6516SJed Brown }
3915edf6516SJed Brown 
39218be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
39318be62a5SSatish Balay {
39418be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
39518be62a5SSatish Balay   PetscErrorCode ierr;
39618be62a5SSatish Balay 
39718be62a5SSatish Balay   PetscFunctionBegin;
39828f357e3SAlex Fikl   ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr);
39918be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
4008fe8eb27SJed Brown   if (shell->dshift) {
4018fe8eb27SJed Brown     ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr);
4028fe8eb27SJed Brown   } else {
40318be62a5SSatish Balay     ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
40418be62a5SSatish Balay   }
4058fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
4068fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
40718be62a5SSatish Balay   PetscFunctionReturn(0);
40818be62a5SSatish Balay }
40918be62a5SSatish Balay 
410f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
411ef66eb69SBarry Smith {
412ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4138fe8eb27SJed Brown   PetscErrorCode ierr;
414b24ad042SBarry Smith 
415ef66eb69SBarry Smith   PetscFunctionBegin;
4168fe8eb27SJed Brown   if (shell->left || shell->right || shell->dshift) {
4178fe8eb27SJed Brown     if (!shell->dshift) {
4188fe8eb27SJed Brown       if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);}
4198fe8eb27SJed Brown       shell->dshift = shell->dshift_owned;
4208fe8eb27SJed Brown       ierr          = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr);
4218fe8eb27SJed Brown     } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
4228fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
4238fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
4248fe8eb27SJed Brown   } else shell->vshift += a;
4258fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
426ef66eb69SBarry Smith   PetscFunctionReturn(0);
427ef66eb69SBarry Smith }
428ef66eb69SBarry Smith 
4296464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
4306464bf51SAlex Fikl {
4316464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
4326464bf51SAlex Fikl   PetscErrorCode ierr;
4336464bf51SAlex Fikl 
4346464bf51SAlex Fikl   PetscFunctionBegin;
4356464bf51SAlex Fikl   if (shell->vscale != 1.0 || shell->left || shell->right) {
4366464bf51SAlex Fikl     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported on scaled matrix");
4376464bf51SAlex Fikl   }
4386464bf51SAlex Fikl 
43928f357e3SAlex Fikl   if (shell->ops->diagonalset) {
44028f357e3SAlex Fikl     ierr = (*shell->ops->diagonalset)(A,D,ins);CHKERRQ(ierr);
44128f357e3SAlex Fikl   }
4426464bf51SAlex Fikl   shell->vshift = 0.0;
4436464bf51SAlex Fikl   if (shell->dshift) { ierr = VecZeroEntries(shell->dshift);CHKERRQ(ierr); }
4446464bf51SAlex Fikl   PetscFunctionReturn(0);
4456464bf51SAlex Fikl }
4466464bf51SAlex Fikl 
447f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
448ef66eb69SBarry Smith {
449ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4508fe8eb27SJed Brown   PetscErrorCode ierr;
451b24ad042SBarry Smith 
452ef66eb69SBarry Smith   PetscFunctionBegin;
453f4df32b1SMatthew Knepley   shell->vscale *= a;
4542205254eSKarl Rupp   if (shell->dshift) {
4552205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
4562205254eSKarl Rupp   } else shell->vshift *= a;
4578fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
4588fe8eb27SJed Brown   PetscFunctionReturn(0);
45918be62a5SSatish Balay }
4608fe8eb27SJed Brown 
4618fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
4628fe8eb27SJed Brown {
4638fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4648fe8eb27SJed Brown   PetscErrorCode ierr;
4658fe8eb27SJed Brown 
4668fe8eb27SJed Brown   PetscFunctionBegin;
4678fe8eb27SJed Brown   if (left) {
4688fe8eb27SJed Brown     if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);}
4692205254eSKarl Rupp     if (shell->left) {
4702205254eSKarl Rupp       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
4712205254eSKarl Rupp     } else {
4728fe8eb27SJed Brown       shell->left = shell->left_owned;
4738fe8eb27SJed Brown       ierr        = VecCopy(left,shell->left);CHKERRQ(ierr);
47418be62a5SSatish Balay     }
475ef66eb69SBarry Smith   }
4768fe8eb27SJed Brown   if (right) {
4778fe8eb27SJed Brown     if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);}
4782205254eSKarl Rupp     if (shell->right) {
4792205254eSKarl Rupp       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
4802205254eSKarl Rupp     } else {
4818fe8eb27SJed Brown       shell->right = shell->right_owned;
4828fe8eb27SJed Brown       ierr         = VecCopy(right,shell->right);CHKERRQ(ierr);
4838fe8eb27SJed Brown     }
4848fe8eb27SJed Brown   }
4858fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
486ef66eb69SBarry Smith   PetscFunctionReturn(0);
487ef66eb69SBarry Smith }
488ef66eb69SBarry Smith 
489dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
490ef66eb69SBarry Smith {
491ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
492ef66eb69SBarry Smith 
493ef66eb69SBarry Smith   PetscFunctionBegin;
4948fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
495ef66eb69SBarry Smith     shell->vshift = 0.0;
496ef66eb69SBarry Smith     shell->vscale = 1.0;
4970298fd71SBarry Smith     shell->dshift = NULL;
4980298fd71SBarry Smith     shell->left   = NULL;
4990298fd71SBarry Smith     shell->right  = NULL;
50028f357e3SAlex Fikl     if (shell->ops->mult) {
50128f357e3SAlex Fikl       Y->ops->mult = shell->ops->mult;
50228f357e3SAlex Fikl       shell->ops->mult  = NULL;
5037fb7d96aSJed Brown     }
50428f357e3SAlex Fikl     if (shell->ops->multtranspose) {
50528f357e3SAlex Fikl       Y->ops->multtranspose = shell->ops->multtranspose;
50628f357e3SAlex Fikl       shell->ops->multtranspose = NULL;
5077fb7d96aSJed Brown     }
50828f357e3SAlex Fikl     if (shell->ops->getdiagonal) {
50928f357e3SAlex Fikl       Y->ops->getdiagonal = shell->ops->getdiagonal;
51028f357e3SAlex Fikl       shell->ops->getdiagonal = NULL;
51128f357e3SAlex Fikl     }
51228f357e3SAlex Fikl     if (shell->ops->copy) {
51328f357e3SAlex Fikl       Y->ops->copy = shell->ops->copy;
51428f357e3SAlex Fikl       shell->ops->copy = NULL;
5157fb7d96aSJed Brown     }
5167fb7d96aSJed Brown     shell->usingscaled = PETSC_FALSE;
517ef66eb69SBarry Smith   }
518ef66eb69SBarry Smith   PetscFunctionReturn(0);
519ef66eb69SBarry Smith }
520ef66eb69SBarry Smith 
521cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
522b951964fSBarry Smith 
5233b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
5243b49f96aSBarry Smith {
5253b49f96aSBarry Smith   PetscFunctionBegin;
5263b49f96aSBarry Smith   *missing = PETSC_FALSE;
5273b49f96aSBarry Smith   PetscFunctionReturn(0);
5283b49f96aSBarry Smith }
5293b49f96aSBarry Smith 
53009dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
53120563c6bSBarry Smith                                        0,
53220563c6bSBarry Smith                                        0,
53320563c6bSBarry Smith                                        0,
53497304618SKris Buschelman                                 /* 4*/ 0,
53520563c6bSBarry Smith                                        0,
536b951964fSBarry Smith                                        0,
537b951964fSBarry Smith                                        0,
538b951964fSBarry Smith                                        0,
539b951964fSBarry Smith                                        0,
54097304618SKris Buschelman                                 /*10*/ 0,
541b951964fSBarry Smith                                        0,
542b951964fSBarry Smith                                        0,
543b951964fSBarry Smith                                        0,
544b951964fSBarry Smith                                        0,
54597304618SKris Buschelman                                 /*15*/ 0,
546b951964fSBarry Smith                                        0,
547b951964fSBarry Smith                                        0,
5488fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
549b951964fSBarry Smith                                        0,
55097304618SKris Buschelman                                 /*20*/ 0,
551ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
552b951964fSBarry Smith                                        0,
553b951964fSBarry Smith                                        0,
554d519adbfSMatthew Knepley                                 /*24*/ 0,
555b951964fSBarry Smith                                        0,
556b951964fSBarry Smith                                        0,
557b951964fSBarry Smith                                        0,
558b951964fSBarry Smith                                        0,
559d519adbfSMatthew Knepley                                 /*29*/ 0,
560b951964fSBarry Smith                                        0,
561273d9f13SBarry Smith                                        0,
562b951964fSBarry Smith                                        0,
563b951964fSBarry Smith                                        0,
564*cb8c8a77SBarry Smith                                 /*34*/ MatDuplicate_Shell,
565b951964fSBarry Smith                                        0,
566b951964fSBarry Smith                                        0,
56709dc0095SBarry Smith                                        0,
56809dc0095SBarry Smith                                        0,
569d519adbfSMatthew Knepley                                 /*39*/ 0,
57009dc0095SBarry Smith                                        0,
57109dc0095SBarry Smith                                        0,
57209dc0095SBarry Smith                                        0,
573*cb8c8a77SBarry Smith                                        MatCopy_Shell,
574d519adbfSMatthew Knepley                                 /*44*/ 0,
575ef66eb69SBarry Smith                                        MatScale_Shell,
576ef66eb69SBarry Smith                                        MatShift_Shell,
5776464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
57809dc0095SBarry Smith                                        0,
579f73d5cc4SBarry Smith                                 /*49*/ 0,
58009dc0095SBarry Smith                                        0,
58109dc0095SBarry Smith                                        0,
58209dc0095SBarry Smith                                        0,
58309dc0095SBarry Smith                                        0,
584d519adbfSMatthew Knepley                                 /*54*/ 0,
58509dc0095SBarry Smith                                        0,
58609dc0095SBarry Smith                                        0,
58709dc0095SBarry Smith                                        0,
58809dc0095SBarry Smith                                        0,
589d519adbfSMatthew Knepley                                 /*59*/ 0,
590b9b97703SBarry Smith                                        MatDestroy_Shell,
59109dc0095SBarry Smith                                        0,
592357abbc8SBarry Smith                                        0,
593273d9f13SBarry Smith                                        0,
594d519adbfSMatthew Knepley                                 /*64*/ 0,
595273d9f13SBarry Smith                                        0,
596273d9f13SBarry Smith                                        0,
597273d9f13SBarry Smith                                        0,
598273d9f13SBarry Smith                                        0,
599d519adbfSMatthew Knepley                                 /*69*/ 0,
600273d9f13SBarry Smith                                        0,
601c87e5d42SMatthew Knepley                                        MatConvert_Shell,
602273d9f13SBarry Smith                                        0,
60397304618SKris Buschelman                                        0,
604d519adbfSMatthew Knepley                                 /*74*/ 0,
60597304618SKris Buschelman                                        0,
60697304618SKris Buschelman                                        0,
60797304618SKris Buschelman                                        0,
60897304618SKris Buschelman                                        0,
609d519adbfSMatthew Knepley                                 /*79*/ 0,
61097304618SKris Buschelman                                        0,
61197304618SKris Buschelman                                        0,
61297304618SKris Buschelman                                        0,
613865e5f61SKris Buschelman                                        0,
614d519adbfSMatthew Knepley                                 /*84*/ 0,
615865e5f61SKris Buschelman                                        0,
616865e5f61SKris Buschelman                                        0,
617865e5f61SKris Buschelman                                        0,
618865e5f61SKris Buschelman                                        0,
619d519adbfSMatthew Knepley                                 /*89*/ 0,
620865e5f61SKris Buschelman                                        0,
621865e5f61SKris Buschelman                                        0,
622865e5f61SKris Buschelman                                        0,
623865e5f61SKris Buschelman                                        0,
624d519adbfSMatthew Knepley                                 /*94*/ 0,
625865e5f61SKris Buschelman                                        0,
626865e5f61SKris Buschelman                                        0,
6273964eb88SJed Brown                                        0,
6283964eb88SJed Brown                                        0,
6293964eb88SJed Brown                                 /*99*/ 0,
6303964eb88SJed Brown                                        0,
6313964eb88SJed Brown                                        0,
6323964eb88SJed Brown                                        0,
6333964eb88SJed Brown                                        0,
6343964eb88SJed Brown                                /*104*/ 0,
6353964eb88SJed Brown                                        0,
6363964eb88SJed Brown                                        0,
6373964eb88SJed Brown                                        0,
6383964eb88SJed Brown                                        0,
6393964eb88SJed Brown                                /*109*/ 0,
6403964eb88SJed Brown                                        0,
6413964eb88SJed Brown                                        0,
6423964eb88SJed Brown                                        0,
6433b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
6443964eb88SJed Brown                                /*114*/ 0,
6453964eb88SJed Brown                                        0,
6463964eb88SJed Brown                                        0,
6473964eb88SJed Brown                                        0,
6483964eb88SJed Brown                                        0,
6493964eb88SJed Brown                                /*119*/ 0,
6503964eb88SJed Brown                                        0,
6513964eb88SJed Brown                                        0,
6523964eb88SJed Brown                                        0,
6533964eb88SJed Brown                                        0,
6543964eb88SJed Brown                                /*124*/ 0,
6553964eb88SJed Brown                                        0,
6563964eb88SJed Brown                                        0,
6573964eb88SJed Brown                                        0,
6583964eb88SJed Brown                                        0,
6593964eb88SJed Brown                                /*129*/ 0,
6603964eb88SJed Brown                                        0,
6613964eb88SJed Brown                                        0,
6623964eb88SJed Brown                                        0,
6633964eb88SJed Brown                                        0,
6643964eb88SJed Brown                                /*134*/ 0,
6653964eb88SJed Brown                                        0,
6663964eb88SJed Brown                                        0,
6673964eb88SJed Brown                                        0,
6683964eb88SJed Brown                                        0,
6693964eb88SJed Brown                                /*139*/ 0,
670f9426fe0SMark Adams                                        0,
6713964eb88SJed Brown                                        0
6723964eb88SJed Brown };
673273d9f13SBarry Smith 
6740bad9183SKris Buschelman /*MC
675fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
6760bad9183SKris Buschelman 
6770bad9183SKris Buschelman   Level: advanced
6780bad9183SKris Buschelman 
6790bad9183SKris Buschelman .seealso: MatCreateShell
6800bad9183SKris Buschelman M*/
6810bad9183SKris Buschelman 
6828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
683273d9f13SBarry Smith {
684273d9f13SBarry Smith   Mat_Shell      *b;
685dfbe8321SBarry Smith   PetscErrorCode ierr;
686273d9f13SBarry Smith 
687273d9f13SBarry Smith   PetscFunctionBegin;
688273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
689273d9f13SBarry Smith 
690b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
691273d9f13SBarry Smith   A->data = (void*)b;
692273d9f13SBarry Smith 
69326283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
69426283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
695273d9f13SBarry Smith 
696273d9f13SBarry Smith   b->ctx           = 0;
697ef66eb69SBarry Smith   b->vshift        = 0.0;
698ef66eb69SBarry Smith   b->vscale        = 1.0;
699273d9f13SBarry Smith   A->assembled     = PETSC_TRUE;
700210f0121SBarry Smith   A->preallocated  = PETSC_FALSE;
7012205254eSKarl Rupp 
70217667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
703273d9f13SBarry Smith   PetscFunctionReturn(0);
704273d9f13SBarry Smith }
705e51e0e81SBarry Smith 
7064b828684SBarry Smith /*@C
707052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
708ff756334SLois Curfman McInnes    private data storage format.
709e51e0e81SBarry Smith 
710c7fcc2eaSBarry Smith   Collective on MPI_Comm
711c7fcc2eaSBarry Smith 
712e51e0e81SBarry Smith    Input Parameters:
713c7fcc2eaSBarry Smith +  comm - MPI communicator
714c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
715c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
716c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
717c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
718c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
719e51e0e81SBarry Smith 
720ff756334SLois Curfman McInnes    Output Parameter:
72144cd7ae7SLois Curfman McInnes .  A - the matrix
722e51e0e81SBarry Smith 
723ff2fd236SBarry Smith    Level: advanced
724ff2fd236SBarry Smith 
725f39d1f56SLois Curfman McInnes   Usage:
7267b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
727f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
728c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
729f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
730f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
731f39d1f56SLois Curfman McInnes 
732ff756334SLois Curfman McInnes    Notes:
733ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
734ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
735ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
736e51e0e81SBarry Smith 
737daf670e6SBarry Smith    Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this
738daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
739daf670e6SBarry Smith     in as the ctx argument.
740f38a8d56SBarry Smith 
741f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
742f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
743645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
744645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
745645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
746645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
747645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
748645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
749645985a0SLois Curfman McInnes    For example,
750f39d1f56SLois Curfman McInnes 
751f39d1f56SLois Curfman McInnes $
752f39d1f56SLois Curfman McInnes $     Vec x, y
7537b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
754645985a0SLois Curfman McInnes $     Mat A
755f39d1f56SLois Curfman McInnes $
756c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
757c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
758f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
759c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
760c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
761c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
762645985a0SLois Curfman McInnes $     MatMult(A,x,y);
763645985a0SLois Curfman McInnes $     MatDestroy(A);
764f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
765645985a0SLois Curfman McInnes $
766e51e0e81SBarry Smith 
7670b627109SLois Curfman McInnes .keywords: matrix, shell, create
7680b627109SLois Curfman McInnes 
769ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
770e51e0e81SBarry Smith @*/
7717087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
772e51e0e81SBarry Smith {
773dfbe8321SBarry Smith   PetscErrorCode ierr;
774ed3cc1f0SBarry Smith 
7753a40ed3dSBarry Smith   PetscFunctionBegin;
776f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
777f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
778273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
779273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
7800fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
781273d9f13SBarry Smith   PetscFunctionReturn(0);
782c7fcc2eaSBarry Smith }
783c7fcc2eaSBarry Smith 
784c6866cfdSSatish Balay /*@
785273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
786c7fcc2eaSBarry Smith 
7873f9fe445SBarry Smith    Logically Collective on Mat
788c7fcc2eaSBarry Smith 
789273d9f13SBarry Smith     Input Parameters:
790273d9f13SBarry Smith +   mat - the shell matrix
791273d9f13SBarry Smith -   ctx - the context
792273d9f13SBarry Smith 
793273d9f13SBarry Smith    Level: advanced
794273d9f13SBarry Smith 
795daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
796daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
797273d9f13SBarry Smith 
798273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
7990bc0a719SBarry Smith @*/
8007087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
801273d9f13SBarry Smith {
802273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
803dfbe8321SBarry Smith   PetscErrorCode ierr;
804ace3abfcSBarry Smith   PetscBool      flg;
805273d9f13SBarry Smith 
806273d9f13SBarry Smith   PetscFunctionBegin;
8070700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
808251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
809273d9f13SBarry Smith   if (flg) {
810273d9f13SBarry Smith     shell->ctx = ctx;
811ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
8123a40ed3dSBarry Smith   PetscFunctionReturn(0);
813e51e0e81SBarry Smith }
814e51e0e81SBarry Smith 
815c16cb8f2SBarry Smith /*@C
8163a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
8173a3eedf2SBarry Smith                            a shell matrix.
818e51e0e81SBarry Smith 
8193f9fe445SBarry Smith    Logically Collective on Mat
820fee21e36SBarry Smith 
821c7fcc2eaSBarry Smith     Input Parameters:
822c7fcc2eaSBarry Smith +   mat - the shell matrix
823c7fcc2eaSBarry Smith .   op - the name of the operation
824c7fcc2eaSBarry Smith -   f - the function that provides the operation.
825c7fcc2eaSBarry Smith 
82615091d37SBarry Smith    Level: advanced
82715091d37SBarry Smith 
828fae171e0SBarry Smith     Usage:
829e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
830f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
831c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
8320b627109SLois Curfman McInnes 
833a62d957aSLois Curfman McInnes     Notes:
834e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
8351c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
836a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
8371c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
838a62d957aSLois Curfman McInnes 
839a4d85fd7SBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
840deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
841deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
842deebb3c3SLois Curfman McInnes     routines, e.g.,
843a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
844a62d957aSLois Curfman McInnes 
8454aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
8464aa34b0aSBarry Smith     nonzero on failure.
8474aa34b0aSBarry Smith 
848a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
849a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
850a62d957aSLois Curfman McInnes     set by MatCreateShell().
851a62d957aSLois Curfman McInnes 
8522a7a6963SBarry Smith     Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
853501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
854501d9185SBarry Smith 
855a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
856a62d957aSLois Curfman McInnes 
857ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
858e51e0e81SBarry Smith @*/
8597087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
860e51e0e81SBarry Smith {
861dfbe8321SBarry Smith   PetscErrorCode ierr;
862ace3abfcSBarry Smith   PetscBool      flg;
863273d9f13SBarry Smith 
8643a40ed3dSBarry Smith   PetscFunctionBegin;
8650700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
8665edf6516SJed Brown   switch (op) {
8675edf6516SJed Brown   case MATOP_DESTROY:
868251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
869273d9f13SBarry Smith     if (flg) {
870a62d957aSLois Curfman McInnes       Mat_Shell *shell = (Mat_Shell*)mat->data;
87128f357e3SAlex Fikl       shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
8726849ba73SBarry Smith     } else mat->ops->destroy = (PetscErrorCode (*)(Mat))f;
8735edf6516SJed Brown     break;
8746464bf51SAlex Fikl   case MATOP_DIAGONAL_SET:
8756464bf51SAlex Fikl     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
8766464bf51SAlex Fikl     if (flg) {
8776464bf51SAlex Fikl       Mat_Shell *shell = (Mat_Shell*)mat->data;
87828f357e3SAlex Fikl       shell->ops->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f;
8796464bf51SAlex Fikl     } else {
8806464bf51SAlex Fikl       mat->ops->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f;
8816464bf51SAlex Fikl     }
8826464bf51SAlex Fikl     break;
8835edf6516SJed Brown   case MATOP_VIEW:
88440a5a6c4SBarry Smith     if (!mat->ops->viewnative) {
88540a5a6c4SBarry Smith       mat->ops->viewnative = mat->ops->view;
88640a5a6c4SBarry Smith     }
8875edf6516SJed Brown     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
8885edf6516SJed Brown     break;
8895edf6516SJed Brown   case MATOP_MULT:
8905edf6516SJed Brown     mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
891875c1aa2SLisandro Dalcin     if (!mat->ops->multadd) {
892875c1aa2SLisandro Dalcin       ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
893875c1aa2SLisandro Dalcin       if (flg) mat->ops->multadd = MatMultAdd_Shell;
894875c1aa2SLisandro Dalcin     }
8955edf6516SJed Brown     break;
8965edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
8975edf6516SJed Brown     mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
898875c1aa2SLisandro Dalcin     if (!mat->ops->multtransposeadd) {
899875c1aa2SLisandro Dalcin       ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
900875c1aa2SLisandro Dalcin       if (flg) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell;
901875c1aa2SLisandro Dalcin     }
9025edf6516SJed Brown     break;
9035edf6516SJed Brown   default:
9045edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
905a62d957aSLois Curfman McInnes   }
9063a40ed3dSBarry Smith   PetscFunctionReturn(0);
907e51e0e81SBarry Smith }
908f0479e8cSBarry Smith 
909d4bb536fSBarry Smith /*@C
910d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
911d4bb536fSBarry Smith 
912c7fcc2eaSBarry Smith     Not Collective
913c7fcc2eaSBarry Smith 
914d4bb536fSBarry Smith     Input Parameters:
915c7fcc2eaSBarry Smith +   mat - the shell matrix
916c7fcc2eaSBarry Smith -   op - the name of the operation
917d4bb536fSBarry Smith 
918d4bb536fSBarry Smith     Output Parameter:
919d4bb536fSBarry Smith .   f - the function that provides the operation.
920d4bb536fSBarry Smith 
92115091d37SBarry Smith     Level: advanced
92215091d37SBarry Smith 
923d4bb536fSBarry Smith     Notes:
924e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
925d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
926d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
927d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
928d4bb536fSBarry Smith 
929d4bb536fSBarry Smith     All user-provided functions have the same calling
930d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
931d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
932d4bb536fSBarry Smith     routines, e.g.,
933d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
934d4bb536fSBarry Smith 
935d4bb536fSBarry Smith     Within each user-defined routine, the user should call
936d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
937d4bb536fSBarry Smith     set by MatCreateShell().
938d4bb536fSBarry Smith 
939d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
940d4bb536fSBarry Smith 
941ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
942d4bb536fSBarry Smith @*/
9437087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
944d4bb536fSBarry Smith {
945dfbe8321SBarry Smith   PetscErrorCode ierr;
946ace3abfcSBarry Smith   PetscBool      flg;
947273d9f13SBarry Smith 
9483a40ed3dSBarry Smith   PetscFunctionBegin;
9490700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
95028f357e3SAlex Fikl   switch (op) {
95128f357e3SAlex Fikl   case MATOP_DESTROY:
952251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
953273d9f13SBarry Smith     if (flg) {
954d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
95528f357e3SAlex Fikl       *f = (void (*)(void))shell->ops->destroy;
956c7fcc2eaSBarry Smith     } else {
957c134de8dSSatish Balay       *f = (void (*)(void))mat->ops->destroy;
958d4bb536fSBarry Smith     }
95928f357e3SAlex Fikl     break;
96028f357e3SAlex Fikl   case MATOP_DIAGONAL_SET:
96128f357e3SAlex Fikl     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
96228f357e3SAlex Fikl     if (flg) {
96328f357e3SAlex Fikl       Mat_Shell *shell = (Mat_Shell*)mat->data;
96428f357e3SAlex Fikl       *f = (void (*)(void))shell->ops->diagonalset;
965c7fcc2eaSBarry Smith     } else {
96628f357e3SAlex Fikl       *f = (void (*)(void))mat->ops->diagonalset;
96728f357e3SAlex Fikl     }
96828f357e3SAlex Fikl     break;
96928f357e3SAlex Fikl   case MATOP_VIEW:
97028f357e3SAlex Fikl     *f = (void (*)(void))mat->ops->view;
97128f357e3SAlex Fikl     break;
97228f357e3SAlex Fikl   default:
973c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
974d4bb536fSBarry Smith   }
9753a40ed3dSBarry Smith   PetscFunctionReturn(0);
976d4bb536fSBarry Smith }
977d4bb536fSBarry Smith 
978