xref: /petsc/src/mat/impls/shell/shell.c (revision 28f357e33c484785ab2e6ff28544298727ab203a)
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 
10*28f357e3SAlex Fikl struct _MatShellOps {
11*28f357e3SAlex Fikl   /*   0 */
126849ba73SBarry Smith   PetscErrorCode (*mult)(Mat,Vec,Vec);
13*28f357e3SAlex Fikl   /*   5 */
1418be62a5SSatish Balay   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
15*28f357e3SAlex Fikl   /*  10 */
16*28f357e3SAlex Fikl   /*  15 */
1718be62a5SSatish Balay   PetscErrorCode (*getdiagonal)(Mat,Vec);
18*28f357e3SAlex Fikl   /*  20 */
19*28f357e3SAlex Fikl   /*  24 */
20*28f357e3SAlex Fikl   /*  29 */
21*28f357e3SAlex Fikl   /*  34 */
22*28f357e3SAlex Fikl   /*  39 */
23*28f357e3SAlex Fikl   PetscErrorCode (*copy)(Mat,Mat,MatStructure);
24*28f357e3SAlex Fikl   /*  44 */
256464bf51SAlex Fikl   PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode);
26*28f357e3SAlex Fikl   /*  49 */
27*28f357e3SAlex Fikl   /*  54 */
28*28f357e3SAlex Fikl   /*  59 */
29*28f357e3SAlex Fikl   PetscErrorCode (*destroy)(Mat);
30*28f357e3SAlex Fikl   /*  64 */
31*28f357e3SAlex Fikl   /*  69 */
32*28f357e3SAlex Fikl   /*  74 */
33*28f357e3SAlex Fikl   /*  79 */
34*28f357e3SAlex Fikl   /*  84 */
35*28f357e3SAlex Fikl   /*  89 */
36*28f357e3SAlex Fikl   /*  94 */
37*28f357e3SAlex Fikl   /*  99 */
38*28f357e3SAlex Fikl   /* 104 */
39*28f357e3SAlex Fikl   /* 109 */
40*28f357e3SAlex Fikl   /* 114 */
41*28f357e3SAlex Fikl   /* 119 */
42*28f357e3SAlex Fikl   /* 124 */
43*28f357e3SAlex Fikl   /* 129 */
44*28f357e3SAlex Fikl   /* 134 */
45*28f357e3SAlex Fikl   /* 139 */
46*28f357e3SAlex Fikl   /* 144 */
47*28f357e3SAlex Fikl };
48*28f357e3SAlex Fikl 
49*28f357e3SAlex Fikl typedef struct {
50*28f357e3SAlex 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);
95*28f357e3SAlex Fikl static PetscErrorCode MatCopy_Shell(Mat,Mat,MatStructure);
968fe8eb27SJed Brown 
978fe8eb27SJed Brown #undef __FUNCT__
988fe8eb27SJed Brown #define __FUNCT__ "MatShellUseScaledMethods"
998fe8eb27SJed Brown static PetscErrorCode MatShellUseScaledMethods(Mat Y)
1008fe8eb27SJed Brown {
1018fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)Y->data;
1028fe8eb27SJed Brown 
1038fe8eb27SJed Brown   PetscFunctionBegin;
1048fe8eb27SJed Brown   if (shell->usingscaled) PetscFunctionReturn(0);
105*28f357e3SAlex Fikl   shell->ops->mult  = Y->ops->mult;
1068fe8eb27SJed Brown   Y->ops->mult = MatMult_Shell;
1078fe8eb27SJed Brown   if (Y->ops->multtranspose) {
108*28f357e3SAlex Fikl     shell->ops->multtranspose  = Y->ops->multtranspose;
1098fe8eb27SJed Brown     Y->ops->multtranspose = MatMultTranspose_Shell;
1108fe8eb27SJed Brown   }
1118fe8eb27SJed Brown   if (Y->ops->getdiagonal) {
112*28f357e3SAlex Fikl     shell->ops->getdiagonal  = Y->ops->getdiagonal;
1138fe8eb27SJed Brown     Y->ops->getdiagonal = MatGetDiagonal_Shell;
1148fe8eb27SJed Brown   }
115*28f357e3SAlex Fikl   if (Y->ops->copy) {
116*28f357e3SAlex Fikl     shell->ops->copy  = Y->ops->copy;
117*28f357e3SAlex Fikl     Y->ops->copy = MatCopy_Shell;
118*28f357e3SAlex Fikl   }
1198fe8eb27SJed Brown   shell->usingscaled = PETSC_TRUE;
1208fe8eb27SJed Brown   PetscFunctionReturn(0);
1218fe8eb27SJed Brown }
1228fe8eb27SJed Brown 
1238fe8eb27SJed Brown #undef __FUNCT__
1248fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleLeft"
1258fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
1268fe8eb27SJed Brown {
1278fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1288fe8eb27SJed Brown   PetscErrorCode ierr;
1298fe8eb27SJed Brown 
1308fe8eb27SJed Brown   PetscFunctionBegin;
1310298fd71SBarry Smith   *xx = NULL;
1328fe8eb27SJed Brown   if (!shell->left) {
1338fe8eb27SJed Brown     *xx = x;
1348fe8eb27SJed Brown   } else {
1358fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
1368fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
1378fe8eb27SJed Brown     *xx  = shell->left_work;
1388fe8eb27SJed Brown   }
1398fe8eb27SJed Brown   PetscFunctionReturn(0);
1408fe8eb27SJed Brown }
1418fe8eb27SJed Brown 
1428fe8eb27SJed Brown #undef __FUNCT__
1438fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleRight"
1448fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
1458fe8eb27SJed Brown {
1468fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1478fe8eb27SJed Brown   PetscErrorCode ierr;
1488fe8eb27SJed Brown 
1498fe8eb27SJed Brown   PetscFunctionBegin;
1500298fd71SBarry Smith   *xx = NULL;
1518fe8eb27SJed Brown   if (!shell->right) {
1528fe8eb27SJed Brown     *xx = x;
1538fe8eb27SJed Brown   } else {
1548fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
1558fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
1568fe8eb27SJed Brown     *xx  = shell->right_work;
1578fe8eb27SJed Brown   }
1588fe8eb27SJed Brown   PetscFunctionReturn(0);
1598fe8eb27SJed Brown }
1608fe8eb27SJed Brown 
1618fe8eb27SJed Brown #undef __FUNCT__
1628fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleLeft"
1638fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
1648fe8eb27SJed Brown {
1658fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1668fe8eb27SJed Brown   PetscErrorCode ierr;
1678fe8eb27SJed Brown 
1688fe8eb27SJed Brown   PetscFunctionBegin;
1698fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
1708fe8eb27SJed Brown   PetscFunctionReturn(0);
1718fe8eb27SJed Brown }
1728fe8eb27SJed Brown 
1738fe8eb27SJed Brown #undef __FUNCT__
1748fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleRight"
1758fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
1768fe8eb27SJed Brown {
1778fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1788fe8eb27SJed Brown   PetscErrorCode ierr;
1798fe8eb27SJed Brown 
1808fe8eb27SJed Brown   PetscFunctionBegin;
1818fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
1828fe8eb27SJed Brown   PetscFunctionReturn(0);
1838fe8eb27SJed Brown }
1848fe8eb27SJed Brown 
1858fe8eb27SJed Brown #undef __FUNCT__
1868fe8eb27SJed Brown #define __FUNCT__ "MatShellShiftAndScale"
1878fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
1888fe8eb27SJed Brown {
1898fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1908fe8eb27SJed Brown   PetscErrorCode ierr;
1918fe8eb27SJed Brown 
1928fe8eb27SJed Brown   PetscFunctionBegin;
1938fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
1948fe8eb27SJed Brown     PetscInt          i,m;
1958fe8eb27SJed Brown     const PetscScalar *x,*d;
1968fe8eb27SJed Brown     PetscScalar       *y;
1978fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
1988fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1998fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
2008fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
2018fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
2028fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
2038fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
2048fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
205880538faSHong Zhang   } else if (PetscAbsScalar(shell->vshift) != 0) {
2068fe8eb27SJed Brown     ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr);
207027c5fb5SJed Brown   } else if (shell->vscale != 1.0) {
208027c5fb5SJed Brown     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
2098fe8eb27SJed Brown   }
2108fe8eb27SJed Brown   PetscFunctionReturn(0);
2118fe8eb27SJed Brown }
212e51e0e81SBarry Smith 
213711e205bSSatish Balay #undef __FUNCT__
214711e205bSSatish Balay #define __FUNCT__ "MatShellGetContext"
2159d225801SJed Brown /*@
216a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
217b4fd4287SBarry Smith 
218c7fcc2eaSBarry Smith     Not Collective
219c7fcc2eaSBarry Smith 
220b4fd4287SBarry Smith     Input Parameter:
221b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
222b4fd4287SBarry Smith 
223b4fd4287SBarry Smith     Output Parameter:
224b4fd4287SBarry Smith .   ctx - the user provided context
225b4fd4287SBarry Smith 
22615091d37SBarry Smith     Level: advanced
22715091d37SBarry Smith 
228daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
229daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
230a62d957aSLois Curfman McInnes 
231a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context
232a62d957aSLois Curfman McInnes 
233ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
2340bc0a719SBarry Smith @*/
2358fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
236b4fd4287SBarry Smith {
237dfbe8321SBarry Smith   PetscErrorCode ierr;
238ace3abfcSBarry Smith   PetscBool      flg;
239273d9f13SBarry Smith 
2403a40ed3dSBarry Smith   PetscFunctionBegin;
2410700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2424482741eSBarry Smith   PetscValidPointer(ctx,2);
243251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
244940183f3SBarry Smith   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
245ce94432eSBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
2463a40ed3dSBarry Smith   PetscFunctionReturn(0);
247b4fd4287SBarry Smith }
248b4fd4287SBarry Smith 
249711e205bSSatish Balay #undef __FUNCT__
250711e205bSSatish Balay #define __FUNCT__ "MatDestroy_Shell"
251dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
252e51e0e81SBarry Smith {
253dfbe8321SBarry Smith   PetscErrorCode ierr;
254bf0cc555SLisandro Dalcin   Mat_Shell      *shell = (Mat_Shell*)mat->data;
255ed3cc1f0SBarry Smith 
2563a40ed3dSBarry Smith   PetscFunctionBegin;
257*28f357e3SAlex Fikl   if (shell->ops->destroy) {
258*28f357e3SAlex Fikl     ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr);
259bf0cc555SLisandro Dalcin   }
2608fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr);
2618fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr);
2628fe8eb27SJed Brown   ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr);
2638fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
2648fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
2655edf6516SJed Brown   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
2665edf6516SJed Brown   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
267bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
2683a40ed3dSBarry Smith   PetscFunctionReturn(0);
269e51e0e81SBarry Smith }
270e51e0e81SBarry Smith 
271711e205bSSatish Balay #undef __FUNCT__
2727fabda3fSAlex Fikl #define __FUNCT__ "MatCopy_Shell"
2737fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
2747fabda3fSAlex Fikl {
275*28f357e3SAlex Fikl   Mat_Shell       *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
276*28f357e3SAlex Fikl   PetscBool       matflg,shellflg;
2777fabda3fSAlex Fikl   PetscErrorCode  ierr;
2787fabda3fSAlex Fikl 
2797fabda3fSAlex Fikl   PetscFunctionBegin;
280*28f357e3SAlex Fikl   ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr);
281*28f357e3SAlex Fikl   if(!matflg) { SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL"); }
2827fabda3fSAlex Fikl 
2837fabda3fSAlex Fikl   ierr = MatShellUseScaledMethods(B);CHKERRQ(ierr);
284*28f357e3SAlex Fikl   ierr = PetscMemcmp(A->ops,B->ops,sizeof(struct _MatOps),&matflg);CHKERRQ(ierr);
285*28f357e3SAlex Fikl   ierr = PetscMemcmp(shellA->ops,shellB->ops,sizeof(struct _MatShellOps),&shellflg);CHKERRQ(ierr);
286*28f357e3SAlex Fikl   if (!matflg || !shellflg) { SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Matrices cannot be copied because they have different operations"); }
2877fabda3fSAlex Fikl 
288*28f357e3SAlex Fikl   ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr);
2897fabda3fSAlex Fikl   shellB->vscale = shellA->vscale;
2907fabda3fSAlex Fikl   shellB->vshift = shellA->vshift;
2917fabda3fSAlex Fikl   if (shellA->dshift_owned) {
2927fabda3fSAlex Fikl     if (!shellB->dshift_owned) {
2937fabda3fSAlex Fikl       ierr = VecDuplicate(shellA->dshift_owned,&shellB->dshift_owned);CHKERRQ(ierr);
2947fabda3fSAlex Fikl     }
2957fabda3fSAlex Fikl     ierr = VecCopy(shellA->dshift_owned,shellB->dshift_owned);CHKERRQ(ierr);
2967fabda3fSAlex Fikl     shellB->dshift = shellB->dshift_owned;
2977fabda3fSAlex Fikl   } else {
2987fabda3fSAlex Fikl     shellB->dshift = NULL;
2997fabda3fSAlex Fikl   }
3007fabda3fSAlex Fikl   if (shellA->left_owned) {
3017fabda3fSAlex Fikl     if (!shellB->left_owned) {
3027fabda3fSAlex Fikl       ierr = VecDuplicate(shellA->left_owned,&shellB->left_owned);CHKERRQ(ierr);
3037fabda3fSAlex Fikl     }
3047fabda3fSAlex Fikl     ierr = VecCopy(shellA->left_owned,shellB->left_owned);CHKERRQ(ierr);
3057fabda3fSAlex Fikl     shellB->left = shellB->left_owned;
3067fabda3fSAlex Fikl   } else {
3077fabda3fSAlex Fikl     shellB->left = NULL;
3087fabda3fSAlex Fikl   }
3097fabda3fSAlex Fikl   if (shellA->right_owned) {
3107fabda3fSAlex Fikl     if (!shellB->right_owned) {
3117fabda3fSAlex Fikl       ierr = VecDuplicate(shellA->right_owned,&shellB->right_owned);CHKERRQ(ierr);
3127fabda3fSAlex Fikl     }
3137fabda3fSAlex Fikl     ierr = VecCopy(shellA->right_owned,shellB->right_owned);CHKERRQ(ierr);
3147fabda3fSAlex Fikl     shellB->right = shellB->right_owned;
3157fabda3fSAlex Fikl   } else {
3167fabda3fSAlex Fikl     shellB->right = NULL;
3177fabda3fSAlex Fikl   }
3187fabda3fSAlex Fikl   PetscFunctionReturn(0);
3197fabda3fSAlex Fikl }
3207fabda3fSAlex Fikl 
3217fabda3fSAlex Fikl #undef __FUNCT__
322ef66eb69SBarry Smith #define __FUNCT__ "MatMult_Shell"
323dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
324ef66eb69SBarry Smith {
325ef66eb69SBarry Smith   Mat_Shell        *shell = (Mat_Shell*)A->data;
326dfbe8321SBarry Smith   PetscErrorCode   ierr;
32725578ef6SJed Brown   Vec              xx;
328e3f487b0SBarry Smith   PetscObjectState instate,outstate;
329ef66eb69SBarry Smith 
330ef66eb69SBarry Smith   PetscFunctionBegin;
3318fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
332e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
333*28f357e3SAlex Fikl   ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr);
334e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
335e3f487b0SBarry Smith   if (instate == outstate) {
336e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
337e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
338e3f487b0SBarry Smith   }
3398fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
3408fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
341ef66eb69SBarry Smith   PetscFunctionReturn(0);
342ef66eb69SBarry Smith }
343ef66eb69SBarry Smith 
344ef66eb69SBarry Smith #undef __FUNCT__
3455edf6516SJed Brown #define __FUNCT__ "MatMultAdd_Shell"
3465edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3475edf6516SJed Brown {
3485edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3495edf6516SJed Brown   PetscErrorCode ierr;
3505edf6516SJed Brown 
3515edf6516SJed Brown   PetscFunctionBegin;
3525edf6516SJed Brown   if (y == z) {
3535edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
3545edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
355b8430d49SBarry Smith     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
3565edf6516SJed Brown   } else {
3575edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
3585edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3595edf6516SJed Brown   }
3605edf6516SJed Brown   PetscFunctionReturn(0);
3615edf6516SJed Brown }
3625edf6516SJed Brown 
3635edf6516SJed Brown #undef __FUNCT__
36418be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell"
36518be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
36618be62a5SSatish Balay {
36718be62a5SSatish Balay   Mat_Shell        *shell = (Mat_Shell*)A->data;
36818be62a5SSatish Balay   PetscErrorCode   ierr;
36925578ef6SJed Brown   Vec              xx;
370e3f487b0SBarry Smith   PetscObjectState instate,outstate;
37118be62a5SSatish Balay 
37218be62a5SSatish Balay   PetscFunctionBegin;
3738fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
374e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
375*28f357e3SAlex Fikl   ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr);
376e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
377e3f487b0SBarry Smith   if (instate == outstate) {
378e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
379e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
380e3f487b0SBarry Smith   }
3818fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
3828fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
38318be62a5SSatish Balay   PetscFunctionReturn(0);
38418be62a5SSatish Balay }
38518be62a5SSatish Balay 
38618be62a5SSatish Balay #undef __FUNCT__
3875edf6516SJed Brown #define __FUNCT__ "MatMultTransposeAdd_Shell"
3885edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3895edf6516SJed Brown {
3905edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3915edf6516SJed Brown   PetscErrorCode ierr;
3925edf6516SJed Brown 
3935edf6516SJed Brown   PetscFunctionBegin;
3945edf6516SJed Brown   if (y == z) {
3955edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
3965edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
3975edf6516SJed Brown     ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr);
3985edf6516SJed Brown   } else {
3995edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
4005edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
4015edf6516SJed Brown   }
4025edf6516SJed Brown   PetscFunctionReturn(0);
4035edf6516SJed Brown }
4045edf6516SJed Brown 
4055edf6516SJed Brown #undef __FUNCT__
40618be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell"
40718be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
40818be62a5SSatish Balay {
40918be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
41018be62a5SSatish Balay   PetscErrorCode ierr;
41118be62a5SSatish Balay 
41218be62a5SSatish Balay   PetscFunctionBegin;
413*28f357e3SAlex Fikl   ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr);
41418be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
4158fe8eb27SJed Brown   if (shell->dshift) {
4168fe8eb27SJed Brown     ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr);
4178fe8eb27SJed Brown   } else {
41818be62a5SSatish Balay     ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
41918be62a5SSatish Balay   }
4208fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
4218fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
42218be62a5SSatish Balay   PetscFunctionReturn(0);
42318be62a5SSatish Balay }
42418be62a5SSatish Balay 
42518be62a5SSatish Balay #undef __FUNCT__
426ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell"
427f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
428ef66eb69SBarry Smith {
429ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4308fe8eb27SJed Brown   PetscErrorCode ierr;
431b24ad042SBarry Smith 
432ef66eb69SBarry Smith   PetscFunctionBegin;
4338fe8eb27SJed Brown   if (shell->left || shell->right || shell->dshift) {
4348fe8eb27SJed Brown     if (!shell->dshift) {
4358fe8eb27SJed Brown       if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);}
4368fe8eb27SJed Brown       shell->dshift = shell->dshift_owned;
4378fe8eb27SJed Brown       ierr          = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr);
4388fe8eb27SJed Brown     } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
4398fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
4408fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
4418fe8eb27SJed Brown   } else shell->vshift += a;
4428fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
443ef66eb69SBarry Smith   PetscFunctionReturn(0);
444ef66eb69SBarry Smith }
445ef66eb69SBarry Smith 
446ef66eb69SBarry Smith #undef __FUNCT__
4476464bf51SAlex Fikl #define __FUNCT__ "MatDiagonalSet_Shell"
4486464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
4496464bf51SAlex Fikl {
4506464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
4516464bf51SAlex Fikl   PetscErrorCode ierr;
4526464bf51SAlex Fikl 
4536464bf51SAlex Fikl   PetscFunctionBegin;
4546464bf51SAlex Fikl   if (shell->vscale != 1.0 || shell->left || shell->right) {
4556464bf51SAlex Fikl     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported on scaled matrix");
4566464bf51SAlex Fikl   }
4576464bf51SAlex Fikl 
458*28f357e3SAlex Fikl   if (shell->ops->diagonalset) {
459*28f357e3SAlex Fikl     ierr = (*shell->ops->diagonalset)(A,D,ins);CHKERRQ(ierr);
460*28f357e3SAlex Fikl   }
4616464bf51SAlex Fikl   shell->vshift = 0.0;
4626464bf51SAlex Fikl   if (shell->dshift) { ierr = VecZeroEntries(shell->dshift);CHKERRQ(ierr); }
4636464bf51SAlex Fikl   PetscFunctionReturn(0);
4646464bf51SAlex Fikl }
4656464bf51SAlex Fikl 
4666464bf51SAlex Fikl #undef __FUNCT__
467ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell"
468f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
469ef66eb69SBarry Smith {
470ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4718fe8eb27SJed Brown   PetscErrorCode ierr;
472b24ad042SBarry Smith 
473ef66eb69SBarry Smith   PetscFunctionBegin;
474f4df32b1SMatthew Knepley   shell->vscale *= a;
4752205254eSKarl Rupp   if (shell->dshift) {
4762205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
4772205254eSKarl Rupp   } else shell->vshift *= a;
4788fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
4798fe8eb27SJed Brown   PetscFunctionReturn(0);
48018be62a5SSatish Balay }
4818fe8eb27SJed Brown 
4828fe8eb27SJed Brown #undef __FUNCT__
4838fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell"
4848fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
4858fe8eb27SJed Brown {
4868fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4878fe8eb27SJed Brown   PetscErrorCode ierr;
4888fe8eb27SJed Brown 
4898fe8eb27SJed Brown   PetscFunctionBegin;
4908fe8eb27SJed Brown   if (left) {
4918fe8eb27SJed Brown     if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);}
4922205254eSKarl Rupp     if (shell->left) {
4932205254eSKarl Rupp       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
4942205254eSKarl Rupp     } else {
4958fe8eb27SJed Brown       shell->left = shell->left_owned;
4968fe8eb27SJed Brown       ierr        = VecCopy(left,shell->left);CHKERRQ(ierr);
49718be62a5SSatish Balay     }
498ef66eb69SBarry Smith   }
4998fe8eb27SJed Brown   if (right) {
5008fe8eb27SJed Brown     if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);}
5012205254eSKarl Rupp     if (shell->right) {
5022205254eSKarl Rupp       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
5032205254eSKarl Rupp     } else {
5048fe8eb27SJed Brown       shell->right = shell->right_owned;
5058fe8eb27SJed Brown       ierr         = VecCopy(right,shell->right);CHKERRQ(ierr);
5068fe8eb27SJed Brown     }
5078fe8eb27SJed Brown   }
5088fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
509ef66eb69SBarry Smith   PetscFunctionReturn(0);
510ef66eb69SBarry Smith }
511ef66eb69SBarry Smith 
512ef66eb69SBarry Smith #undef __FUNCT__
513ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell"
514dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
515ef66eb69SBarry Smith {
516ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
517ef66eb69SBarry Smith 
518ef66eb69SBarry Smith   PetscFunctionBegin;
5198fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
520ef66eb69SBarry Smith     shell->vshift = 0.0;
521ef66eb69SBarry Smith     shell->vscale = 1.0;
5220298fd71SBarry Smith     shell->dshift = NULL;
5230298fd71SBarry Smith     shell->left   = NULL;
5240298fd71SBarry Smith     shell->right  = NULL;
525*28f357e3SAlex Fikl     if (shell->ops->mult) {
526*28f357e3SAlex Fikl       Y->ops->mult = shell->ops->mult;
527*28f357e3SAlex Fikl       shell->ops->mult  = NULL;
5287fb7d96aSJed Brown     }
529*28f357e3SAlex Fikl     if (shell->ops->multtranspose) {
530*28f357e3SAlex Fikl       Y->ops->multtranspose = shell->ops->multtranspose;
531*28f357e3SAlex Fikl       shell->ops->multtranspose = NULL;
5327fb7d96aSJed Brown     }
533*28f357e3SAlex Fikl     if (shell->ops->getdiagonal) {
534*28f357e3SAlex Fikl       Y->ops->getdiagonal = shell->ops->getdiagonal;
535*28f357e3SAlex Fikl       shell->ops->getdiagonal = NULL;
536*28f357e3SAlex Fikl     }
537*28f357e3SAlex Fikl     if (shell->ops->copy) {
538*28f357e3SAlex Fikl       Y->ops->copy = shell->ops->copy;
539*28f357e3SAlex Fikl       shell->ops->copy = NULL;
5407fb7d96aSJed Brown     }
5417fb7d96aSJed Brown     shell->usingscaled = PETSC_FALSE;
542ef66eb69SBarry Smith   }
543ef66eb69SBarry Smith   PetscFunctionReturn(0);
544ef66eb69SBarry Smith }
545ef66eb69SBarry Smith 
546cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
547b951964fSBarry Smith 
5483b49f96aSBarry Smith #undef __FUNCT__
5493b49f96aSBarry Smith #define __FUNCT__ "MatMissingDiagonal_Shell"
5503b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
5513b49f96aSBarry Smith {
5523b49f96aSBarry Smith   PetscFunctionBegin;
5533b49f96aSBarry Smith   *missing = PETSC_FALSE;
5543b49f96aSBarry Smith   PetscFunctionReturn(0);
5553b49f96aSBarry Smith }
5563b49f96aSBarry Smith 
55709dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
55820563c6bSBarry Smith                                        0,
55920563c6bSBarry Smith                                        0,
56020563c6bSBarry Smith                                        0,
56197304618SKris Buschelman                                 /* 4*/ 0,
56220563c6bSBarry Smith                                        0,
563b951964fSBarry Smith                                        0,
564b951964fSBarry Smith                                        0,
565b951964fSBarry Smith                                        0,
566b951964fSBarry Smith                                        0,
56797304618SKris Buschelman                                 /*10*/ 0,
568b951964fSBarry Smith                                        0,
569b951964fSBarry Smith                                        0,
570b951964fSBarry Smith                                        0,
571b951964fSBarry Smith                                        0,
57297304618SKris Buschelman                                 /*15*/ 0,
573b951964fSBarry Smith                                        0,
574b951964fSBarry Smith                                        0,
5758fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
576b951964fSBarry Smith                                        0,
57797304618SKris Buschelman                                 /*20*/ 0,
578ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
579b951964fSBarry Smith                                        0,
580b951964fSBarry Smith                                        0,
581d519adbfSMatthew Knepley                                 /*24*/ 0,
582b951964fSBarry Smith                                        0,
583b951964fSBarry Smith                                        0,
584b951964fSBarry Smith                                        0,
585b951964fSBarry Smith                                        0,
586d519adbfSMatthew Knepley                                 /*29*/ 0,
587b951964fSBarry Smith                                        0,
588273d9f13SBarry Smith                                        0,
589b951964fSBarry Smith                                        0,
590b951964fSBarry Smith                                        0,
591d519adbfSMatthew Knepley                                 /*34*/ 0,
592b951964fSBarry Smith                                        0,
593b951964fSBarry Smith                                        0,
59409dc0095SBarry Smith                                        0,
59509dc0095SBarry Smith                                        0,
596d519adbfSMatthew Knepley                                 /*39*/ 0,
59709dc0095SBarry Smith                                        0,
59809dc0095SBarry Smith                                        0,
59909dc0095SBarry Smith                                        0,
600*28f357e3SAlex Fikl                                        0,
601d519adbfSMatthew Knepley                                 /*44*/ 0,
602ef66eb69SBarry Smith                                        MatScale_Shell,
603ef66eb69SBarry Smith                                        MatShift_Shell,
6046464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
60509dc0095SBarry Smith                                        0,
606f73d5cc4SBarry Smith                                 /*49*/ 0,
60709dc0095SBarry Smith                                        0,
60809dc0095SBarry Smith                                        0,
60909dc0095SBarry Smith                                        0,
61009dc0095SBarry Smith                                        0,
611d519adbfSMatthew Knepley                                 /*54*/ 0,
61209dc0095SBarry Smith                                        0,
61309dc0095SBarry Smith                                        0,
61409dc0095SBarry Smith                                        0,
61509dc0095SBarry Smith                                        0,
616d519adbfSMatthew Knepley                                 /*59*/ 0,
617b9b97703SBarry Smith                                        MatDestroy_Shell,
61809dc0095SBarry Smith                                        0,
619357abbc8SBarry Smith                                        0,
620273d9f13SBarry Smith                                        0,
621d519adbfSMatthew Knepley                                 /*64*/ 0,
622273d9f13SBarry Smith                                        0,
623273d9f13SBarry Smith                                        0,
624273d9f13SBarry Smith                                        0,
625273d9f13SBarry Smith                                        0,
626d519adbfSMatthew Knepley                                 /*69*/ 0,
627273d9f13SBarry Smith                                        0,
628c87e5d42SMatthew Knepley                                        MatConvert_Shell,
629273d9f13SBarry Smith                                        0,
63097304618SKris Buschelman                                        0,
631d519adbfSMatthew Knepley                                 /*74*/ 0,
63297304618SKris Buschelman                                        0,
63397304618SKris Buschelman                                        0,
63497304618SKris Buschelman                                        0,
63597304618SKris Buschelman                                        0,
636d519adbfSMatthew Knepley                                 /*79*/ 0,
63797304618SKris Buschelman                                        0,
63897304618SKris Buschelman                                        0,
63997304618SKris Buschelman                                        0,
640865e5f61SKris Buschelman                                        0,
641d519adbfSMatthew Knepley                                 /*84*/ 0,
642865e5f61SKris Buschelman                                        0,
643865e5f61SKris Buschelman                                        0,
644865e5f61SKris Buschelman                                        0,
645865e5f61SKris Buschelman                                        0,
646d519adbfSMatthew Knepley                                 /*89*/ 0,
647865e5f61SKris Buschelman                                        0,
648865e5f61SKris Buschelman                                        0,
649865e5f61SKris Buschelman                                        0,
650865e5f61SKris Buschelman                                        0,
651d519adbfSMatthew Knepley                                 /*94*/ 0,
652865e5f61SKris Buschelman                                        0,
653865e5f61SKris Buschelman                                        0,
6543964eb88SJed Brown                                        0,
6553964eb88SJed Brown                                        0,
6563964eb88SJed Brown                                 /*99*/ 0,
6573964eb88SJed Brown                                        0,
6583964eb88SJed Brown                                        0,
6593964eb88SJed Brown                                        0,
6603964eb88SJed Brown                                        0,
6613964eb88SJed Brown                                /*104*/ 0,
6623964eb88SJed Brown                                        0,
6633964eb88SJed Brown                                        0,
6643964eb88SJed Brown                                        0,
6653964eb88SJed Brown                                        0,
6663964eb88SJed Brown                                /*109*/ 0,
6673964eb88SJed Brown                                        0,
6683964eb88SJed Brown                                        0,
6693964eb88SJed Brown                                        0,
6703b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
6713964eb88SJed Brown                                /*114*/ 0,
6723964eb88SJed Brown                                        0,
6733964eb88SJed Brown                                        0,
6743964eb88SJed Brown                                        0,
6753964eb88SJed Brown                                        0,
6763964eb88SJed Brown                                /*119*/ 0,
6773964eb88SJed Brown                                        0,
6783964eb88SJed Brown                                        0,
6793964eb88SJed Brown                                        0,
6803964eb88SJed Brown                                        0,
6813964eb88SJed Brown                                /*124*/ 0,
6823964eb88SJed Brown                                        0,
6833964eb88SJed Brown                                        0,
6843964eb88SJed Brown                                        0,
6853964eb88SJed Brown                                        0,
6863964eb88SJed Brown                                /*129*/ 0,
6873964eb88SJed Brown                                        0,
6883964eb88SJed Brown                                        0,
6893964eb88SJed Brown                                        0,
6903964eb88SJed Brown                                        0,
6913964eb88SJed Brown                                /*134*/ 0,
6923964eb88SJed Brown                                        0,
6933964eb88SJed Brown                                        0,
6943964eb88SJed Brown                                        0,
6953964eb88SJed Brown                                        0,
6963964eb88SJed Brown                                /*139*/ 0,
697f9426fe0SMark Adams                                        0,
6983964eb88SJed Brown                                        0
6993964eb88SJed Brown };
700273d9f13SBarry Smith 
7010bad9183SKris Buschelman /*MC
702fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
7030bad9183SKris Buschelman 
7040bad9183SKris Buschelman   Level: advanced
7050bad9183SKris Buschelman 
7060bad9183SKris Buschelman .seealso: MatCreateShell
7070bad9183SKris Buschelman M*/
7080bad9183SKris Buschelman 
709711e205bSSatish Balay #undef __FUNCT__
710711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell"
7118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
712273d9f13SBarry Smith {
713273d9f13SBarry Smith   Mat_Shell      *b;
714dfbe8321SBarry Smith   PetscErrorCode ierr;
715273d9f13SBarry Smith 
716273d9f13SBarry Smith   PetscFunctionBegin;
717273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
718273d9f13SBarry Smith 
719b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
720273d9f13SBarry Smith   A->data = (void*)b;
721273d9f13SBarry Smith 
72226283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
72326283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
724273d9f13SBarry Smith 
725273d9f13SBarry Smith   b->ctx           = 0;
726ef66eb69SBarry Smith   b->vshift        = 0.0;
727ef66eb69SBarry Smith   b->vscale        = 1.0;
728273d9f13SBarry Smith   A->assembled     = PETSC_TRUE;
729210f0121SBarry Smith   A->preallocated  = PETSC_FALSE;
7302205254eSKarl Rupp 
73117667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
732273d9f13SBarry Smith   PetscFunctionReturn(0);
733273d9f13SBarry Smith }
734e51e0e81SBarry Smith 
735711e205bSSatish Balay #undef __FUNCT__
736711e205bSSatish Balay #define __FUNCT__ "MatCreateShell"
7374b828684SBarry Smith /*@C
738052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
739ff756334SLois Curfman McInnes    private data storage format.
740e51e0e81SBarry Smith 
741c7fcc2eaSBarry Smith   Collective on MPI_Comm
742c7fcc2eaSBarry Smith 
743e51e0e81SBarry Smith    Input Parameters:
744c7fcc2eaSBarry Smith +  comm - MPI communicator
745c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
746c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
747c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
748c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
749c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
750e51e0e81SBarry Smith 
751ff756334SLois Curfman McInnes    Output Parameter:
75244cd7ae7SLois Curfman McInnes .  A - the matrix
753e51e0e81SBarry Smith 
754ff2fd236SBarry Smith    Level: advanced
755ff2fd236SBarry Smith 
756f39d1f56SLois Curfman McInnes   Usage:
7577b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
758f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
759c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
760f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
761f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
762f39d1f56SLois Curfman McInnes 
763ff756334SLois Curfman McInnes    Notes:
764ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
765ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
766ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
767e51e0e81SBarry Smith 
768daf670e6SBarry Smith    Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this
769daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
770daf670e6SBarry Smith     in as the ctx argument.
771f38a8d56SBarry Smith 
772f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
773f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
774645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
775645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
776645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
777645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
778645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
779645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
780645985a0SLois Curfman McInnes    For example,
781f39d1f56SLois Curfman McInnes 
782f39d1f56SLois Curfman McInnes $
783f39d1f56SLois Curfman McInnes $     Vec x, y
7847b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
785645985a0SLois Curfman McInnes $     Mat A
786f39d1f56SLois Curfman McInnes $
787c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
788c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
789f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
790c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
791c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
792c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
793645985a0SLois Curfman McInnes $     MatMult(A,x,y);
794645985a0SLois Curfman McInnes $     MatDestroy(A);
795f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
796645985a0SLois Curfman McInnes $
797e51e0e81SBarry Smith 
7980b627109SLois Curfman McInnes .keywords: matrix, shell, create
7990b627109SLois Curfman McInnes 
800ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
801e51e0e81SBarry Smith @*/
8027087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
803e51e0e81SBarry Smith {
804dfbe8321SBarry Smith   PetscErrorCode ierr;
805ed3cc1f0SBarry Smith 
8063a40ed3dSBarry Smith   PetscFunctionBegin;
807f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
808f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
809273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
810273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
8110fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
812273d9f13SBarry Smith   PetscFunctionReturn(0);
813c7fcc2eaSBarry Smith }
814c7fcc2eaSBarry Smith 
815711e205bSSatish Balay #undef __FUNCT__
816711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext"
817c6866cfdSSatish Balay /*@
818273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
819c7fcc2eaSBarry Smith 
8203f9fe445SBarry Smith    Logically Collective on Mat
821c7fcc2eaSBarry Smith 
822273d9f13SBarry Smith     Input Parameters:
823273d9f13SBarry Smith +   mat - the shell matrix
824273d9f13SBarry Smith -   ctx - the context
825273d9f13SBarry Smith 
826273d9f13SBarry Smith    Level: advanced
827273d9f13SBarry Smith 
828daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
829daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
830273d9f13SBarry Smith 
831273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
8320bc0a719SBarry Smith @*/
8337087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
834273d9f13SBarry Smith {
835273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
836dfbe8321SBarry Smith   PetscErrorCode ierr;
837ace3abfcSBarry Smith   PetscBool      flg;
838273d9f13SBarry Smith 
839273d9f13SBarry Smith   PetscFunctionBegin;
8400700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
841251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
842273d9f13SBarry Smith   if (flg) {
843273d9f13SBarry Smith     shell->ctx = ctx;
844ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
8453a40ed3dSBarry Smith   PetscFunctionReturn(0);
846e51e0e81SBarry Smith }
847e51e0e81SBarry Smith 
848711e205bSSatish Balay #undef __FUNCT__
849711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation"
850c16cb8f2SBarry Smith /*@C
8513a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
8523a3eedf2SBarry Smith                            a shell matrix.
853e51e0e81SBarry Smith 
8543f9fe445SBarry Smith    Logically Collective on Mat
855fee21e36SBarry Smith 
856c7fcc2eaSBarry Smith     Input Parameters:
857c7fcc2eaSBarry Smith +   mat - the shell matrix
858c7fcc2eaSBarry Smith .   op - the name of the operation
859c7fcc2eaSBarry Smith -   f - the function that provides the operation.
860c7fcc2eaSBarry Smith 
86115091d37SBarry Smith    Level: advanced
86215091d37SBarry Smith 
863fae171e0SBarry Smith     Usage:
864e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
865f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
866c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
8670b627109SLois Curfman McInnes 
868a62d957aSLois Curfman McInnes     Notes:
869e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
8701c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
871a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
8721c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
873a62d957aSLois Curfman McInnes 
87425296bd5SBarry Smith     All user-provided functions (execept for MATOP_DESTROY) should have the same calling
875deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
876deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
877deebb3c3SLois Curfman McInnes     routines, e.g.,
878a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
879a62d957aSLois Curfman McInnes 
8804aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
8814aa34b0aSBarry Smith     nonzero on failure.
8824aa34b0aSBarry Smith 
883a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
884a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
885a62d957aSLois Curfman McInnes     set by MatCreateShell().
886a62d957aSLois Curfman McInnes 
8872a7a6963SBarry Smith     Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
888501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
889501d9185SBarry Smith 
890a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
891a62d957aSLois Curfman McInnes 
892ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
893e51e0e81SBarry Smith @*/
8947087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
895e51e0e81SBarry Smith {
896dfbe8321SBarry Smith   PetscErrorCode ierr;
897ace3abfcSBarry Smith   PetscBool      flg;
898273d9f13SBarry Smith 
8993a40ed3dSBarry Smith   PetscFunctionBegin;
9000700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9015edf6516SJed Brown   switch (op) {
9025edf6516SJed Brown   case MATOP_DESTROY:
903251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
904273d9f13SBarry Smith     if (flg) {
905a62d957aSLois Curfman McInnes       Mat_Shell *shell = (Mat_Shell*)mat->data;
906*28f357e3SAlex Fikl       shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
9076849ba73SBarry Smith     } else mat->ops->destroy = (PetscErrorCode (*)(Mat))f;
9085edf6516SJed Brown     break;
9096464bf51SAlex Fikl   case MATOP_DIAGONAL_SET:
9106464bf51SAlex Fikl     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
9116464bf51SAlex Fikl     if (flg) {
9126464bf51SAlex Fikl       Mat_Shell *shell = (Mat_Shell*)mat->data;
913*28f357e3SAlex Fikl       shell->ops->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f;
9146464bf51SAlex Fikl     } else {
9156464bf51SAlex Fikl       mat->ops->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f;
9166464bf51SAlex Fikl     }
9176464bf51SAlex Fikl     break;
9185edf6516SJed Brown   case MATOP_VIEW:
9195edf6516SJed Brown     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
9205edf6516SJed Brown     break;
9215edf6516SJed Brown   case MATOP_MULT:
9225edf6516SJed Brown     mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
9235edf6516SJed Brown     if (!mat->ops->multadd) mat->ops->multadd = MatMultAdd_Shell;
9245edf6516SJed Brown     break;
9255edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
9265edf6516SJed Brown     mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
9275edf6516SJed Brown     if (!mat->ops->multtransposeadd) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell;
9285edf6516SJed Brown     break;
9295edf6516SJed Brown   default:
9305edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
931a62d957aSLois Curfman McInnes   }
9323a40ed3dSBarry Smith   PetscFunctionReturn(0);
933e51e0e81SBarry Smith }
934f0479e8cSBarry Smith 
935711e205bSSatish Balay #undef __FUNCT__
936711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation"
937d4bb536fSBarry Smith /*@C
938d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
939d4bb536fSBarry Smith 
940c7fcc2eaSBarry Smith     Not Collective
941c7fcc2eaSBarry Smith 
942d4bb536fSBarry Smith     Input Parameters:
943c7fcc2eaSBarry Smith +   mat - the shell matrix
944c7fcc2eaSBarry Smith -   op - the name of the operation
945d4bb536fSBarry Smith 
946d4bb536fSBarry Smith     Output Parameter:
947d4bb536fSBarry Smith .   f - the function that provides the operation.
948d4bb536fSBarry Smith 
94915091d37SBarry Smith     Level: advanced
95015091d37SBarry Smith 
951d4bb536fSBarry Smith     Notes:
952e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
953d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
954d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
955d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
956d4bb536fSBarry Smith 
957d4bb536fSBarry Smith     All user-provided functions have the same calling
958d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
959d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
960d4bb536fSBarry Smith     routines, e.g.,
961d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
962d4bb536fSBarry Smith 
963d4bb536fSBarry Smith     Within each user-defined routine, the user should call
964d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
965d4bb536fSBarry Smith     set by MatCreateShell().
966d4bb536fSBarry Smith 
967d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
968d4bb536fSBarry Smith 
969ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
970d4bb536fSBarry Smith @*/
9717087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
972d4bb536fSBarry Smith {
973dfbe8321SBarry Smith   PetscErrorCode ierr;
974ace3abfcSBarry Smith   PetscBool      flg;
975273d9f13SBarry Smith 
9763a40ed3dSBarry Smith   PetscFunctionBegin;
9770700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
978*28f357e3SAlex Fikl   switch (op) {
979*28f357e3SAlex Fikl   case MATOP_DESTROY:
980251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
981273d9f13SBarry Smith     if (flg) {
982d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
983*28f357e3SAlex Fikl       *f = (void (*)(void))shell->ops->destroy;
984c7fcc2eaSBarry Smith     } else {
985c134de8dSSatish Balay       *f = (void (*)(void))mat->ops->destroy;
986d4bb536fSBarry Smith     }
987*28f357e3SAlex Fikl     break;
988*28f357e3SAlex Fikl   case MATOP_DIAGONAL_SET:
989*28f357e3SAlex Fikl     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
990*28f357e3SAlex Fikl     if (flg) {
991*28f357e3SAlex Fikl       Mat_Shell *shell = (Mat_Shell*)mat->data;
992*28f357e3SAlex Fikl       *f = (void (*)(void))shell->ops->diagonalset;
993c7fcc2eaSBarry Smith     } else {
994*28f357e3SAlex Fikl       *f = (void (*)(void))mat->ops->diagonalset;
995*28f357e3SAlex Fikl     }
996*28f357e3SAlex Fikl     break;
997*28f357e3SAlex Fikl   case MATOP_VIEW:
998*28f357e3SAlex Fikl     *f = (void (*)(void))mat->ops->view;
999*28f357e3SAlex Fikl     break;
1000*28f357e3SAlex Fikl   default:
1001c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
1002d4bb536fSBarry Smith   }
10033a40ed3dSBarry Smith   PetscFunctionReturn(0);
1004d4bb536fSBarry Smith }
1005d4bb536fSBarry Smith 
1006