xref: /petsc/src/mat/impls/shell/shell.c (revision 5259c5a4b12ef8e698ce5dd89d5bafabdb8bd5aa)
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 */
239f137db4SBarry Smith   PetscErrorCode (*axpy)(Mat,PetscScalar,Mat,MatStructure);
2428f357e3SAlex Fikl   PetscErrorCode (*copy)(Mat,Mat,MatStructure);
2528f357e3SAlex Fikl   /*  44 */
266464bf51SAlex Fikl   PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode);
2728f357e3SAlex Fikl   /*  49 */
2828f357e3SAlex Fikl   /*  54 */
2928f357e3SAlex Fikl   /*  59 */
3028f357e3SAlex Fikl   PetscErrorCode (*destroy)(Mat);
3128f357e3SAlex Fikl   /*  64 */
3228f357e3SAlex Fikl   /*  69 */
3328f357e3SAlex Fikl   /*  74 */
3428f357e3SAlex Fikl   /*  79 */
3528f357e3SAlex Fikl   /*  84 */
3628f357e3SAlex Fikl   /*  89 */
3728f357e3SAlex Fikl   /*  94 */
3828f357e3SAlex Fikl   /*  99 */
3928f357e3SAlex Fikl   /* 104 */
4028f357e3SAlex Fikl   /* 109 */
4128f357e3SAlex Fikl   /* 114 */
4228f357e3SAlex Fikl   /* 119 */
4328f357e3SAlex Fikl   /* 124 */
4428f357e3SAlex Fikl   /* 129 */
4528f357e3SAlex Fikl   /* 134 */
4628f357e3SAlex Fikl   /* 139 */
4728f357e3SAlex Fikl   /* 144 */
4828f357e3SAlex Fikl };
4928f357e3SAlex Fikl 
5028f357e3SAlex Fikl typedef struct {
5128f357e3SAlex Fikl   struct _MatShellOps ops[1];
522205254eSKarl Rupp 
53ef66eb69SBarry Smith   PetscScalar vscale,vshift;
548fe8eb27SJed Brown   Vec         dshift;
558fe8eb27SJed Brown   Vec         left,right;
568fe8eb27SJed Brown   Vec         left_work,right_work;
575edf6516SJed Brown   Vec         left_add_work,right_add_work;
589f137db4SBarry Smith   Mat         axpy;
599f137db4SBarry Smith   PetscScalar axpy_vscale;
600c0fd78eSBarry Smith   PetscBool   managescalingshifts;                   /* The user will manage the scaling and shifts for the MATSHELL, not the default */
6120563c6bSBarry Smith   void        *ctx;
6288cf3e7dSBarry Smith } Mat_Shell;
630c0fd78eSBarry Smith 
648fe8eb27SJed Brown /*
650c0fd78eSBarry Smith       xx = diag(left)*x
668fe8eb27SJed Brown */
678fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
688fe8eb27SJed Brown {
698fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
708fe8eb27SJed Brown   PetscErrorCode ierr;
718fe8eb27SJed Brown 
728fe8eb27SJed Brown   PetscFunctionBegin;
730298fd71SBarry Smith   *xx = NULL;
748fe8eb27SJed Brown   if (!shell->left) {
758fe8eb27SJed Brown     *xx = x;
768fe8eb27SJed Brown   } else {
778fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
788fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
798fe8eb27SJed Brown     *xx  = shell->left_work;
808fe8eb27SJed Brown   }
818fe8eb27SJed Brown   PetscFunctionReturn(0);
828fe8eb27SJed Brown }
838fe8eb27SJed Brown 
840c0fd78eSBarry Smith /*
850c0fd78eSBarry Smith      xx = diag(right)*x
860c0fd78eSBarry Smith */
878fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
888fe8eb27SJed Brown {
898fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
908fe8eb27SJed Brown   PetscErrorCode ierr;
918fe8eb27SJed Brown 
928fe8eb27SJed Brown   PetscFunctionBegin;
930298fd71SBarry Smith   *xx = NULL;
948fe8eb27SJed Brown   if (!shell->right) {
958fe8eb27SJed Brown     *xx = x;
968fe8eb27SJed Brown   } else {
978fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
988fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
998fe8eb27SJed Brown     *xx  = shell->right_work;
1008fe8eb27SJed Brown   }
1018fe8eb27SJed Brown   PetscFunctionReturn(0);
1028fe8eb27SJed Brown }
1038fe8eb27SJed Brown 
1040c0fd78eSBarry Smith /*
1050c0fd78eSBarry Smith     x = diag(left)*x
1060c0fd78eSBarry Smith */
1078fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
1088fe8eb27SJed Brown {
1098fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1108fe8eb27SJed Brown   PetscErrorCode ierr;
1118fe8eb27SJed Brown 
1128fe8eb27SJed Brown   PetscFunctionBegin;
1138fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
1148fe8eb27SJed Brown   PetscFunctionReturn(0);
1158fe8eb27SJed Brown }
1168fe8eb27SJed Brown 
1170c0fd78eSBarry Smith /*
1180c0fd78eSBarry Smith     x = diag(right)*x
1190c0fd78eSBarry Smith */
1208fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
1218fe8eb27SJed Brown {
1228fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1238fe8eb27SJed Brown   PetscErrorCode ierr;
1248fe8eb27SJed Brown 
1258fe8eb27SJed Brown   PetscFunctionBegin;
1268fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
1278fe8eb27SJed Brown   PetscFunctionReturn(0);
1288fe8eb27SJed Brown }
1298fe8eb27SJed Brown 
1300c0fd78eSBarry Smith /*
1310c0fd78eSBarry Smith          Y = vscale*Y + diag(dshift)*X + vshift*X
1320c0fd78eSBarry Smith 
1330c0fd78eSBarry Smith          On input Y already contains A*x
1340c0fd78eSBarry Smith */
1358fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
1368fe8eb27SJed Brown {
1378fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1388fe8eb27SJed Brown   PetscErrorCode ierr;
1398fe8eb27SJed Brown 
1408fe8eb27SJed Brown   PetscFunctionBegin;
1418fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
1428fe8eb27SJed Brown     PetscInt          i,m;
1438fe8eb27SJed Brown     const PetscScalar *x,*d;
1448fe8eb27SJed Brown     PetscScalar       *y;
1458fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
1468fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1478fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1488fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1498fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
1508fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1518fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1528fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
1530c0fd78eSBarry Smith   } else {
154027c5fb5SJed Brown     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
1558fe8eb27SJed Brown   }
1560c0fd78eSBarry Smith   if (shell->vshift) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */
1578fe8eb27SJed Brown   PetscFunctionReturn(0);
1588fe8eb27SJed Brown }
159e51e0e81SBarry Smith 
1609d225801SJed Brown /*@
161a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
162b4fd4287SBarry Smith 
163c7fcc2eaSBarry Smith     Not Collective
164c7fcc2eaSBarry Smith 
165b4fd4287SBarry Smith     Input Parameter:
166b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
167b4fd4287SBarry Smith 
168b4fd4287SBarry Smith     Output Parameter:
169b4fd4287SBarry Smith .   ctx - the user provided context
170b4fd4287SBarry Smith 
17115091d37SBarry Smith     Level: advanced
17215091d37SBarry Smith 
173daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
174daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
175a62d957aSLois Curfman McInnes 
176a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context
177a62d957aSLois Curfman McInnes 
178ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
1790bc0a719SBarry Smith @*/
1808fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
181b4fd4287SBarry Smith {
182dfbe8321SBarry Smith   PetscErrorCode ierr;
183ace3abfcSBarry Smith   PetscBool      flg;
184273d9f13SBarry Smith 
1853a40ed3dSBarry Smith   PetscFunctionBegin;
1860700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1874482741eSBarry Smith   PetscValidPointer(ctx,2);
188251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
189940183f3SBarry Smith   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
190ce94432eSBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
1913a40ed3dSBarry Smith   PetscFunctionReturn(0);
192b4fd4287SBarry Smith }
193b4fd4287SBarry Smith 
194dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
195e51e0e81SBarry Smith {
196dfbe8321SBarry Smith   PetscErrorCode ierr;
197bf0cc555SLisandro Dalcin   Mat_Shell      *shell = (Mat_Shell*)mat->data;
198ed3cc1f0SBarry Smith 
1993a40ed3dSBarry Smith   PetscFunctionBegin;
20028f357e3SAlex Fikl   if (shell->ops->destroy) {
20128f357e3SAlex Fikl     ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr);
202bf0cc555SLisandro Dalcin   }
2030c0fd78eSBarry Smith   ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
2040c0fd78eSBarry Smith   ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
2050c0fd78eSBarry Smith   ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
2068fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
2078fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
2085edf6516SJed Brown   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
2095edf6516SJed Brown   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
2109f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
211bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
2123a40ed3dSBarry Smith   PetscFunctionReturn(0);
213e51e0e81SBarry Smith }
214e51e0e81SBarry Smith 
2157fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
2167fabda3fSAlex Fikl {
21728f357e3SAlex Fikl   Mat_Shell       *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
2187fabda3fSAlex Fikl   PetscErrorCode  ierr;
219cb8c8a77SBarry Smith   PetscBool       matflg;
2207fabda3fSAlex Fikl 
2217fabda3fSAlex Fikl   PetscFunctionBegin;
22228f357e3SAlex Fikl   ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr);
223cb8c8a77SBarry Smith   if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL");
2247fabda3fSAlex Fikl 
225cb8c8a77SBarry Smith   ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
226cb8c8a77SBarry Smith   ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
2277fabda3fSAlex Fikl 
228cb8c8a77SBarry Smith   if (shellA->ops->copy) {
22928f357e3SAlex Fikl     ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr);
230cb8c8a77SBarry Smith   }
2317fabda3fSAlex Fikl   shellB->vscale = shellA->vscale;
2327fabda3fSAlex Fikl   shellB->vshift = shellA->vshift;
2330c0fd78eSBarry Smith   if (shellA->dshift) {
2340c0fd78eSBarry Smith     if (!shellB->dshift) {
2350c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr);
2367fabda3fSAlex Fikl     }
2370c0fd78eSBarry Smith     ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr);
2387fabda3fSAlex Fikl   } else {
2399f137db4SBarry Smith     ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr);
2407fabda3fSAlex Fikl   }
2410c0fd78eSBarry Smith   if (shellA->left) {
2420c0fd78eSBarry Smith     if (!shellB->left) {
2430c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr);
2447fabda3fSAlex Fikl     }
2450c0fd78eSBarry Smith     ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr);
2467fabda3fSAlex Fikl   } else {
2479f137db4SBarry Smith     ierr = VecDestroy(&shellB->left);CHKERRQ(ierr);
2487fabda3fSAlex Fikl   }
2490c0fd78eSBarry Smith   if (shellA->right) {
2500c0fd78eSBarry Smith     if (!shellB->right) {
2510c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr);
2527fabda3fSAlex Fikl     }
2530c0fd78eSBarry Smith     ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr);
2547fabda3fSAlex Fikl   } else {
2559f137db4SBarry Smith     ierr = VecDestroy(&shellB->right);CHKERRQ(ierr);
2567fabda3fSAlex Fikl   }
2577fabda3fSAlex Fikl   PetscFunctionReturn(0);
2587fabda3fSAlex Fikl }
2597fabda3fSAlex Fikl 
260cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M)
261cb8c8a77SBarry Smith {
262cb8c8a77SBarry Smith   PetscErrorCode ierr;
263cb8c8a77SBarry Smith   void           *ctx;
264cb8c8a77SBarry Smith 
265cb8c8a77SBarry Smith   PetscFunctionBegin;
266cb8c8a77SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
267cb8c8a77SBarry Smith   ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr);
268cb8c8a77SBarry Smith   ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
269cb8c8a77SBarry Smith   PetscFunctionReturn(0);
270cb8c8a77SBarry Smith }
271cb8c8a77SBarry Smith 
272dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
273ef66eb69SBarry Smith {
274ef66eb69SBarry Smith   Mat_Shell        *shell = (Mat_Shell*)A->data;
275dfbe8321SBarry Smith   PetscErrorCode   ierr;
27625578ef6SJed Brown   Vec              xx;
277e3f487b0SBarry Smith   PetscObjectState instate,outstate;
278ef66eb69SBarry Smith 
279ef66eb69SBarry Smith   PetscFunctionBegin;
2808fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
281e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
28228f357e3SAlex Fikl   ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr);
283e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
284e3f487b0SBarry Smith   if (instate == outstate) {
285e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
286e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
287e3f487b0SBarry Smith   }
2888fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2898fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
2909f137db4SBarry Smith 
2919f137db4SBarry Smith   if (shell->axpy) {
2929f137db4SBarry Smith     if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);}
2939f137db4SBarry Smith     ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr);
2949f137db4SBarry Smith     ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
2959f137db4SBarry Smith   }
296ef66eb69SBarry Smith   PetscFunctionReturn(0);
297ef66eb69SBarry Smith }
298ef66eb69SBarry Smith 
2995edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3005edf6516SJed Brown {
3015edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3025edf6516SJed Brown   PetscErrorCode ierr;
3035edf6516SJed Brown 
3045edf6516SJed Brown   PetscFunctionBegin;
3055edf6516SJed Brown   if (y == z) {
3065edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
3075edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
308b8430d49SBarry Smith     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
3095edf6516SJed Brown   } else {
3105edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
3115edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3125edf6516SJed Brown   }
3135edf6516SJed Brown   PetscFunctionReturn(0);
3145edf6516SJed Brown }
3155edf6516SJed Brown 
31618be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
31718be62a5SSatish Balay {
31818be62a5SSatish Balay   Mat_Shell        *shell = (Mat_Shell*)A->data;
31918be62a5SSatish Balay   PetscErrorCode   ierr;
32025578ef6SJed Brown   Vec              xx;
321e3f487b0SBarry Smith   PetscObjectState instate,outstate;
32218be62a5SSatish Balay 
32318be62a5SSatish Balay   PetscFunctionBegin;
3248fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
325e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
3260c0fd78eSBarry Smith   if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL");
32728f357e3SAlex Fikl   ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr);
328e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
329e3f487b0SBarry Smith   if (instate == outstate) {
330e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
331e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
332e3f487b0SBarry Smith   }
3338fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
3348fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
33518be62a5SSatish Balay   PetscFunctionReturn(0);
33618be62a5SSatish Balay }
33718be62a5SSatish Balay 
3385edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3395edf6516SJed Brown {
3405edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3415edf6516SJed Brown   PetscErrorCode ierr;
3425edf6516SJed Brown 
3435edf6516SJed Brown   PetscFunctionBegin;
3445edf6516SJed Brown   if (y == z) {
3455edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
3465edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
3475edf6516SJed Brown     ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr);
3485edf6516SJed Brown   } else {
3495edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
3505edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3515edf6516SJed Brown   }
3525edf6516SJed Brown   PetscFunctionReturn(0);
3535edf6516SJed Brown }
3545edf6516SJed Brown 
3550c0fd78eSBarry Smith /*
3560c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
3570c0fd78eSBarry Smith */
35818be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
35918be62a5SSatish Balay {
36018be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
36118be62a5SSatish Balay   PetscErrorCode ierr;
36218be62a5SSatish Balay 
36318be62a5SSatish Balay   PetscFunctionBegin;
3640c0fd78eSBarry Smith   if (shell->ops->getdiagonal) {
36528f357e3SAlex Fikl     ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr);
3660c0fd78eSBarry Smith   } else {
3670c0fd78eSBarry Smith     ierr = VecSet(v,0.0);CHKERRQ(ierr);
3680c0fd78eSBarry Smith   }
36918be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
3708fe8eb27SJed Brown   if (shell->dshift) {
3710c0fd78eSBarry Smith     ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr);
37218be62a5SSatish Balay   }
3730c0fd78eSBarry Smith   ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
3748fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
3758fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
37618be62a5SSatish Balay   PetscFunctionReturn(0);
37718be62a5SSatish Balay }
37818be62a5SSatish Balay 
379f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
380ef66eb69SBarry Smith {
381ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3828fe8eb27SJed Brown   PetscErrorCode ierr;
383b24ad042SBarry Smith 
384ef66eb69SBarry Smith   PetscFunctionBegin;
3850c0fd78eSBarry Smith   if (shell->left || shell->right) {
3868fe8eb27SJed Brown     if (!shell->dshift) {
3870c0fd78eSBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
3880c0fd78eSBarry Smith       ierr = VecSet(shell->dshift,a);CHKERRQ(ierr);
3890c0fd78eSBarry Smith     } else {
3900c0fd78eSBarry Smith       if (shell->left)  {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
3910c0fd78eSBarry Smith       if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
3929f137db4SBarry Smith       ierr = VecShift(shell->dshift,a);CHKERRQ(ierr);
3930c0fd78eSBarry Smith     }
3948fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
3958fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
3968fe8eb27SJed Brown   } else shell->vshift += a;
397ef66eb69SBarry Smith   PetscFunctionReturn(0);
398ef66eb69SBarry Smith }
399ef66eb69SBarry Smith 
4006464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
4016464bf51SAlex Fikl {
4026464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
4036464bf51SAlex Fikl   PetscErrorCode ierr;
4046464bf51SAlex Fikl 
4056464bf51SAlex Fikl   PetscFunctionBegin;
4060c0fd78eSBarry Smith   if (ins == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported with INSERT_VALUES");
4070c0fd78eSBarry Smith   if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);}
4080c0fd78eSBarry Smith   if (shell->left || shell->right) {
4099f137db4SBarry Smith     if (!shell->right_work) ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);
4109f137db4SBarry Smith     if (shell->left && shell->right)  {
4119f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
4129f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
4139f137db4SBarry Smith     } else if (shell->left) {
4149f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
4159f137db4SBarry Smith     } else {
4169f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr);
4179f137db4SBarry Smith     }
4189f137db4SBarry Smith     if (!shell->dshift) {
4199f137db4SBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
4209f137db4SBarry Smith       ierr = VecCopy(shell->dshift,shell->right_work);CHKERRQ(ierr);
4219f137db4SBarry Smith     } else {
4229f137db4SBarry Smith       ierr = VecAXPY(shell->dshift,1.0,shell->right_work);CHKERRQ(ierr);
4239f137db4SBarry Smith     }
4240c0fd78eSBarry Smith   } else {
4250c0fd78eSBarry Smith     ierr = VecAXPY(shell->dshift,1.0,D);CHKERRQ(ierr);
4266464bf51SAlex Fikl   }
4276464bf51SAlex Fikl   PetscFunctionReturn(0);
4286464bf51SAlex Fikl }
4296464bf51SAlex Fikl 
430f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
431ef66eb69SBarry Smith {
432ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4338fe8eb27SJed Brown   PetscErrorCode ierr;
434b24ad042SBarry Smith 
435ef66eb69SBarry Smith   PetscFunctionBegin;
436f4df32b1SMatthew Knepley   shell->vscale *= a;
4370c0fd78eSBarry Smith   shell->vshift *= a;
4382205254eSKarl Rupp   if (shell->dshift) {
4392205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
4400c0fd78eSBarry Smith   }
4418fe8eb27SJed Brown   PetscFunctionReturn(0);
44218be62a5SSatish Balay }
4438fe8eb27SJed Brown 
4448fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
4458fe8eb27SJed Brown {
4468fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4478fe8eb27SJed Brown   PetscErrorCode ierr;
4488fe8eb27SJed Brown 
4498fe8eb27SJed Brown   PetscFunctionBegin;
4508fe8eb27SJed Brown   if (left) {
4510c0fd78eSBarry Smith     if (!shell->left) {
4520c0fd78eSBarry Smith       ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr);
4538fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
4540c0fd78eSBarry Smith     } else {
4550c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
45618be62a5SSatish Balay     }
457ef66eb69SBarry Smith   }
4588fe8eb27SJed Brown   if (right) {
4590c0fd78eSBarry Smith     if (!shell->right) {
4600c0fd78eSBarry Smith       ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr);
4618fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
4620c0fd78eSBarry Smith     } else {
4630c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
4648fe8eb27SJed Brown     }
4658fe8eb27SJed Brown   }
466ef66eb69SBarry Smith   PetscFunctionReturn(0);
467ef66eb69SBarry Smith }
468ef66eb69SBarry Smith 
469dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
470ef66eb69SBarry Smith {
471ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4720c0fd78eSBarry Smith   PetscErrorCode ierr;
473ef66eb69SBarry Smith 
474ef66eb69SBarry Smith   PetscFunctionBegin;
4758fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
476ef66eb69SBarry Smith     shell->vshift = 0.0;
477ef66eb69SBarry Smith     shell->vscale = 1.0;
4780c0fd78eSBarry Smith     ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
4790c0fd78eSBarry Smith     ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
4800c0fd78eSBarry Smith     ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
481ef66eb69SBarry Smith   }
482ef66eb69SBarry Smith   PetscFunctionReturn(0);
483ef66eb69SBarry Smith }
484ef66eb69SBarry Smith 
485cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
486b951964fSBarry Smith 
4873b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
4883b49f96aSBarry Smith {
4893b49f96aSBarry Smith   PetscFunctionBegin;
4903b49f96aSBarry Smith   *missing = PETSC_FALSE;
4913b49f96aSBarry Smith   PetscFunctionReturn(0);
4923b49f96aSBarry Smith }
4933b49f96aSBarry Smith 
4949f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
4959f137db4SBarry Smith {
4969f137db4SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4979f137db4SBarry Smith   PetscErrorCode ierr;
4989f137db4SBarry Smith 
4999f137db4SBarry Smith   PetscFunctionBegin;
5009f137db4SBarry Smith   ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr);
5019f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
5029f137db4SBarry Smith   shell->axpy        = X;
5039f137db4SBarry Smith   shell->axpy_vscale = a;
5049f137db4SBarry Smith   PetscFunctionReturn(0);
5059f137db4SBarry Smith }
5069f137db4SBarry Smith 
50709dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
50820563c6bSBarry Smith                                        0,
50920563c6bSBarry Smith                                        0,
5109f137db4SBarry Smith                                        0,
5110c0fd78eSBarry Smith                                 /* 4*/ MatMultAdd_Shell,
5129f137db4SBarry Smith                                        0,
5130c0fd78eSBarry Smith                                        MatMultTransposeAdd_Shell,
514b951964fSBarry Smith                                        0,
515b951964fSBarry Smith                                        0,
516b951964fSBarry Smith                                        0,
51797304618SKris Buschelman                                 /*10*/ 0,
518b951964fSBarry Smith                                        0,
519b951964fSBarry Smith                                        0,
520b951964fSBarry Smith                                        0,
521b951964fSBarry Smith                                        0,
52297304618SKris Buschelman                                 /*15*/ 0,
523b951964fSBarry Smith                                        0,
5240c0fd78eSBarry Smith                                        MatGetDiagonal_Shell,
5258fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
526b951964fSBarry Smith                                        0,
52797304618SKris Buschelman                                 /*20*/ 0,
528ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
529b951964fSBarry Smith                                        0,
530b951964fSBarry Smith                                        0,
531d519adbfSMatthew Knepley                                 /*24*/ 0,
532b951964fSBarry Smith                                        0,
533b951964fSBarry Smith                                        0,
534b951964fSBarry Smith                                        0,
535b951964fSBarry Smith                                        0,
536d519adbfSMatthew Knepley                                 /*29*/ 0,
537b951964fSBarry Smith                                        0,
538273d9f13SBarry Smith                                        0,
539b951964fSBarry Smith                                        0,
540b951964fSBarry Smith                                        0,
541cb8c8a77SBarry Smith                                 /*34*/ MatDuplicate_Shell,
542b951964fSBarry Smith                                        0,
543b951964fSBarry Smith                                        0,
54409dc0095SBarry Smith                                        0,
54509dc0095SBarry Smith                                        0,
5469f137db4SBarry Smith                                 /*39*/ MatAXPY_Shell,
54709dc0095SBarry Smith                                        0,
54809dc0095SBarry Smith                                        0,
54909dc0095SBarry Smith                                        0,
550cb8c8a77SBarry Smith                                        MatCopy_Shell,
551d519adbfSMatthew Knepley                                 /*44*/ 0,
552ef66eb69SBarry Smith                                        MatScale_Shell,
553ef66eb69SBarry Smith                                        MatShift_Shell,
5546464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
55509dc0095SBarry Smith                                        0,
556f73d5cc4SBarry Smith                                 /*49*/ 0,
55709dc0095SBarry Smith                                        0,
55809dc0095SBarry Smith                                        0,
55909dc0095SBarry Smith                                        0,
56009dc0095SBarry Smith                                        0,
561d519adbfSMatthew Knepley                                 /*54*/ 0,
56209dc0095SBarry Smith                                        0,
56309dc0095SBarry Smith                                        0,
56409dc0095SBarry Smith                                        0,
56509dc0095SBarry Smith                                        0,
566d519adbfSMatthew Knepley                                 /*59*/ 0,
567b9b97703SBarry Smith                                        MatDestroy_Shell,
56809dc0095SBarry Smith                                        0,
569357abbc8SBarry Smith                                        0,
570273d9f13SBarry Smith                                        0,
571d519adbfSMatthew Knepley                                 /*64*/ 0,
572273d9f13SBarry Smith                                        0,
573273d9f13SBarry Smith                                        0,
574273d9f13SBarry Smith                                        0,
575273d9f13SBarry Smith                                        0,
576d519adbfSMatthew Knepley                                 /*69*/ 0,
577273d9f13SBarry Smith                                        0,
578c87e5d42SMatthew Knepley                                        MatConvert_Shell,
579273d9f13SBarry Smith                                        0,
58097304618SKris Buschelman                                        0,
581d519adbfSMatthew Knepley                                 /*74*/ 0,
58297304618SKris Buschelman                                        0,
58397304618SKris Buschelman                                        0,
58497304618SKris Buschelman                                        0,
58597304618SKris Buschelman                                        0,
586d519adbfSMatthew Knepley                                 /*79*/ 0,
58797304618SKris Buschelman                                        0,
58897304618SKris Buschelman                                        0,
58997304618SKris Buschelman                                        0,
590865e5f61SKris Buschelman                                        0,
591d519adbfSMatthew Knepley                                 /*84*/ 0,
592865e5f61SKris Buschelman                                        0,
593865e5f61SKris Buschelman                                        0,
594865e5f61SKris Buschelman                                        0,
595865e5f61SKris Buschelman                                        0,
596d519adbfSMatthew Knepley                                 /*89*/ 0,
597865e5f61SKris Buschelman                                        0,
598865e5f61SKris Buschelman                                        0,
599865e5f61SKris Buschelman                                        0,
600865e5f61SKris Buschelman                                        0,
601d519adbfSMatthew Knepley                                 /*94*/ 0,
602865e5f61SKris Buschelman                                        0,
603865e5f61SKris Buschelman                                        0,
6043964eb88SJed Brown                                        0,
6053964eb88SJed Brown                                        0,
6063964eb88SJed Brown                                 /*99*/ 0,
6073964eb88SJed Brown                                        0,
6083964eb88SJed Brown                                        0,
6093964eb88SJed Brown                                        0,
6103964eb88SJed Brown                                        0,
6113964eb88SJed Brown                                /*104*/ 0,
6123964eb88SJed Brown                                        0,
6133964eb88SJed Brown                                        0,
6143964eb88SJed Brown                                        0,
6153964eb88SJed Brown                                        0,
6163964eb88SJed Brown                                /*109*/ 0,
6173964eb88SJed Brown                                        0,
6183964eb88SJed Brown                                        0,
6193964eb88SJed Brown                                        0,
6203b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
6213964eb88SJed Brown                                /*114*/ 0,
6223964eb88SJed Brown                                        0,
6233964eb88SJed Brown                                        0,
6243964eb88SJed Brown                                        0,
6253964eb88SJed Brown                                        0,
6263964eb88SJed Brown                                /*119*/ 0,
6273964eb88SJed Brown                                        0,
6283964eb88SJed Brown                                        0,
6293964eb88SJed Brown                                        0,
6303964eb88SJed Brown                                        0,
6313964eb88SJed Brown                                /*124*/ 0,
6323964eb88SJed Brown                                        0,
6333964eb88SJed Brown                                        0,
6343964eb88SJed Brown                                        0,
6353964eb88SJed Brown                                        0,
6363964eb88SJed Brown                                /*129*/ 0,
6373964eb88SJed Brown                                        0,
6383964eb88SJed Brown                                        0,
6393964eb88SJed Brown                                        0,
6403964eb88SJed Brown                                        0,
6413964eb88SJed Brown                                /*134*/ 0,
6423964eb88SJed Brown                                        0,
6433964eb88SJed Brown                                        0,
6443964eb88SJed Brown                                        0,
6453964eb88SJed Brown                                        0,
6463964eb88SJed Brown                                /*139*/ 0,
647f9426fe0SMark Adams                                        0,
6483964eb88SJed Brown                                        0
6493964eb88SJed Brown };
650273d9f13SBarry Smith 
6510bad9183SKris Buschelman /*MC
652fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
6530bad9183SKris Buschelman 
6540bad9183SKris Buschelman   Level: advanced
6550bad9183SKris Buschelman 
6560c0fd78eSBarry Smith .seealso: MatCreateShell()
6570bad9183SKris Buschelman M*/
6580bad9183SKris Buschelman 
6598cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
660273d9f13SBarry Smith {
661273d9f13SBarry Smith   Mat_Shell      *b;
662dfbe8321SBarry Smith   PetscErrorCode ierr;
663273d9f13SBarry Smith 
664273d9f13SBarry Smith   PetscFunctionBegin;
665273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
666273d9f13SBarry Smith 
667b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
668273d9f13SBarry Smith   A->data = (void*)b;
669273d9f13SBarry Smith 
67026283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
67126283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
672273d9f13SBarry Smith 
673273d9f13SBarry Smith   b->ctx                 = 0;
674ef66eb69SBarry Smith   b->vshift              = 0.0;
675ef66eb69SBarry Smith   b->vscale              = 1.0;
6760c0fd78eSBarry Smith   b->managescalingshifts = PETSC_TRUE;
677273d9f13SBarry Smith   A->assembled           = PETSC_TRUE;
678210f0121SBarry Smith   A->preallocated        = PETSC_FALSE;
6792205254eSKarl Rupp 
68017667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
681273d9f13SBarry Smith   PetscFunctionReturn(0);
682273d9f13SBarry Smith }
683e51e0e81SBarry Smith 
6844b828684SBarry Smith /*@C
685052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
686ff756334SLois Curfman McInnes    private data storage format.
687e51e0e81SBarry Smith 
688c7fcc2eaSBarry Smith   Collective on MPI_Comm
689c7fcc2eaSBarry Smith 
690e51e0e81SBarry Smith    Input Parameters:
691c7fcc2eaSBarry Smith +  comm - MPI communicator
692c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
693c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
694c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
695c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
696c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
697e51e0e81SBarry Smith 
698ff756334SLois Curfman McInnes    Output Parameter:
69944cd7ae7SLois Curfman McInnes .  A - the matrix
700e51e0e81SBarry Smith 
701ff2fd236SBarry Smith    Level: advanced
702ff2fd236SBarry Smith 
703f39d1f56SLois Curfman McInnes   Usage:
7047b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
705f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
706c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
707f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
708f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
709f39d1f56SLois Curfman McInnes 
710ff756334SLois Curfman McInnes    Notes:
711ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
712ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
713ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
714e51e0e81SBarry Smith 
715daf670e6SBarry Smith    Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this
716daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
717daf670e6SBarry Smith     in as the ctx argument.
718f38a8d56SBarry Smith 
719f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
720f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
721645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
722645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
723645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
724645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
725645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
726645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
727645985a0SLois Curfman McInnes    For example,
728f39d1f56SLois Curfman McInnes 
7290c0fd78eSBarry Smith    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), and MatScale() internally so these
7300c0fd78eSBarry Smith    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.
7310c0fd78eSBarry Smith 
732f39d1f56SLois Curfman McInnes $
733f39d1f56SLois Curfman McInnes $     Vec x, y
7347b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
735645985a0SLois Curfman McInnes $     Mat A
736f39d1f56SLois Curfman McInnes $
737c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
738c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
739f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
740c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
741c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
742c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
743645985a0SLois Curfman McInnes $     MatMult(A,x,y);
744645985a0SLois Curfman McInnes $     MatDestroy(A);
745f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
746645985a0SLois Curfman McInnes $
747e51e0e81SBarry Smith 
7480c0fd78eSBarry Smith     For rectangular matrices do all the scalings and shifts make sense?
7490c0fd78eSBarry Smith 
7500c0fd78eSBarry Smith     Developers Notes: Regarding shifting and scaling. The general form is
7510c0fd78eSBarry Smith 
7520c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
7530c0fd78eSBarry Smith 
7540c0fd78eSBarry Smith       The order you apply the operations is important. For example if you have a dshift then
7550c0fd78eSBarry Smith       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
7560c0fd78eSBarry Smith       you get s*vscale*A + diag(shift)
7570c0fd78eSBarry Smith 
7580c0fd78eSBarry Smith           A is the user provided function.
7590c0fd78eSBarry Smith 
7600b627109SLois Curfman McInnes .keywords: matrix, shell, create
7610b627109SLois Curfman McInnes 
7620c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
763e51e0e81SBarry Smith @*/
7647087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
765e51e0e81SBarry Smith {
766dfbe8321SBarry Smith   PetscErrorCode ierr;
767ed3cc1f0SBarry Smith 
7683a40ed3dSBarry Smith   PetscFunctionBegin;
769f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
770f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
771273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
772273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
7730fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
774273d9f13SBarry Smith   PetscFunctionReturn(0);
775c7fcc2eaSBarry Smith }
776c7fcc2eaSBarry Smith 
777c6866cfdSSatish Balay /*@
778273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
779c7fcc2eaSBarry Smith 
7803f9fe445SBarry Smith    Logically Collective on Mat
781c7fcc2eaSBarry Smith 
782273d9f13SBarry Smith     Input Parameters:
783273d9f13SBarry Smith +   mat - the shell matrix
784273d9f13SBarry Smith -   ctx - the context
785273d9f13SBarry Smith 
786273d9f13SBarry Smith    Level: advanced
787273d9f13SBarry Smith 
788daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
789daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
790273d9f13SBarry Smith 
791273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
7920bc0a719SBarry Smith @*/
7937087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
794273d9f13SBarry Smith {
795273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
796dfbe8321SBarry Smith   PetscErrorCode ierr;
797ace3abfcSBarry Smith   PetscBool      flg;
798273d9f13SBarry Smith 
799273d9f13SBarry Smith   PetscFunctionBegin;
8000700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
801251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
802273d9f13SBarry Smith   if (flg) {
803273d9f13SBarry Smith     shell->ctx = ctx;
804ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
8053a40ed3dSBarry Smith   PetscFunctionReturn(0);
806e51e0e81SBarry Smith }
807e51e0e81SBarry Smith 
8080c0fd78eSBarry Smith /*@
8090c0fd78eSBarry Smith     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
8100c0fd78eSBarry Smith           after MatCreateShell()
8110c0fd78eSBarry Smith 
8120c0fd78eSBarry Smith    Logically Collective on Mat
8130c0fd78eSBarry Smith 
8140c0fd78eSBarry Smith     Input Parameter:
8150c0fd78eSBarry Smith .   mat - the shell matrix
8160c0fd78eSBarry Smith 
8170c0fd78eSBarry Smith   Level: advanced
8180c0fd78eSBarry Smith 
8190c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
8200c0fd78eSBarry Smith @*/
8210c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A)
8220c0fd78eSBarry Smith {
8230c0fd78eSBarry Smith   PetscErrorCode ierr;
8240c0fd78eSBarry Smith   Mat_Shell      *shell = (Mat_Shell*)A->data;
8250c0fd78eSBarry Smith   PetscBool      flg;
8260c0fd78eSBarry Smith 
8270c0fd78eSBarry Smith   PetscFunctionBegin;
8280c0fd78eSBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
8290c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr);
8300c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
8310c0fd78eSBarry Smith   shell->managescalingshifts = PETSC_FALSE;
8320c0fd78eSBarry Smith   A->ops->diagonalscale = 0;
8330c0fd78eSBarry Smith   A->ops->scale         = 0;
8340c0fd78eSBarry Smith   A->ops->shift         = 0;
8350c0fd78eSBarry Smith   A->ops->diagonalset   = 0;
8360c0fd78eSBarry Smith   PetscFunctionReturn(0);
8370c0fd78eSBarry Smith }
8380c0fd78eSBarry Smith 
839c16cb8f2SBarry Smith /*@C
840f3b1f45cSBarry Smith     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
841f3b1f45cSBarry Smith 
842f3b1f45cSBarry Smith    Logically Collective on Mat
843f3b1f45cSBarry Smith 
844f3b1f45cSBarry Smith     Input Parameters:
845f3b1f45cSBarry Smith +   mat - the shell matrix
846f3b1f45cSBarry Smith .   f - the function
847f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
848f3b1f45cSBarry Smith -   ctx - an optional context for the function
849f3b1f45cSBarry Smith 
850f3b1f45cSBarry Smith    Output Parameter:
851f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
852f3b1f45cSBarry Smith 
853f3b1f45cSBarry Smith    Options Database:
854f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
855f3b1f45cSBarry Smith 
856f3b1f45cSBarry Smith    Level: advanced
857f3b1f45cSBarry Smith 
858f3b1f45cSBarry Smith    Fortran Notes: Not supported from Fortran
859f3b1f45cSBarry Smith 
860f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
861f3b1f45cSBarry Smith @*/
862f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
863f3b1f45cSBarry Smith {
864f3b1f45cSBarry Smith   PetscErrorCode ierr;
865f3b1f45cSBarry Smith   PetscInt       m,n;
866f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
867f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
86874e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
869f3b1f45cSBarry Smith 
870f3b1f45cSBarry Smith   PetscFunctionBegin;
871f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
872f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr);
873f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
874f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr);
875f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
876f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
877f3b1f45cSBarry Smith 
878f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
879f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mat,&Dmat);CHKERRQ(ierr);
880f3b1f45cSBarry Smith 
881f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
882f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
883f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
884f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
885f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
886f3b1f45cSBarry Smith     flag = PETSC_FALSE;
887f3b1f45cSBarry Smith     if (v) {
888f3b1f45cSBarry Smith       ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm));
889f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
890f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
891f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
892f3b1f45cSBarry Smith     }
893f3b1f45cSBarry Smith   } else if (v) {
894f3b1f45cSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");
895f3b1f45cSBarry Smith   }
896f3b1f45cSBarry Smith   if (flg) *flg = flag;
897f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
898f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
899f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
900f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
901f3b1f45cSBarry Smith   PetscFunctionReturn(0);
902f3b1f45cSBarry Smith }
903f3b1f45cSBarry Smith 
904f3b1f45cSBarry Smith /*@C
905f3b1f45cSBarry Smith     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.
906f3b1f45cSBarry Smith 
907f3b1f45cSBarry Smith    Logically Collective on Mat
908f3b1f45cSBarry Smith 
909f3b1f45cSBarry Smith     Input Parameters:
910f3b1f45cSBarry Smith +   mat - the shell matrix
911f3b1f45cSBarry Smith .   f - the function
912f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
913f3b1f45cSBarry Smith -   ctx - an optional context for the function
914f3b1f45cSBarry Smith 
915f3b1f45cSBarry Smith    Output Parameter:
916f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
917f3b1f45cSBarry Smith 
918f3b1f45cSBarry Smith    Options Database:
919f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
920f3b1f45cSBarry Smith 
921f3b1f45cSBarry Smith    Level: advanced
922f3b1f45cSBarry Smith 
923f3b1f45cSBarry Smith    Fortran Notes: Not supported from Fortran
924f3b1f45cSBarry Smith 
925f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
926f3b1f45cSBarry Smith @*/
927f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
928f3b1f45cSBarry Smith {
929f3b1f45cSBarry Smith   PetscErrorCode ierr;
930f3b1f45cSBarry Smith   Vec            x,y,z;
931f3b1f45cSBarry Smith   PetscInt       m,n,M,N;
932f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
933f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
93474e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
935f3b1f45cSBarry Smith 
936f3b1f45cSBarry Smith   PetscFunctionBegin;
937f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
938f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
939f3b1f45cSBarry Smith   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
940f3b1f45cSBarry Smith   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
941f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
942f3b1f45cSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
943f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
944f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
945f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
946f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
947f3b1f45cSBarry Smith   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
948f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr);
949f3b1f45cSBarry Smith 
950f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
951f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
952f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
953f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
954f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
955f3b1f45cSBarry Smith     flag = PETSC_FALSE;
956f3b1f45cSBarry Smith     if (v) {
957f3b1f45cSBarry Smith       ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm));
958f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
959f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
960f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
961f3b1f45cSBarry Smith     }
962f3b1f45cSBarry Smith   } else if (v) {
963f3b1f45cSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");
964f3b1f45cSBarry Smith   }
965f3b1f45cSBarry Smith   if (flg) *flg = flag;
966f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
967f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
968f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
969f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
970f3b1f45cSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
971f3b1f45cSBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
972f3b1f45cSBarry Smith   ierr = VecDestroy(&z);CHKERRQ(ierr);
973f3b1f45cSBarry Smith   PetscFunctionReturn(0);
974f3b1f45cSBarry Smith }
975f3b1f45cSBarry Smith 
976f3b1f45cSBarry Smith /*@C
9773a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
9783a3eedf2SBarry Smith                            a shell matrix.
9790c0fd78eSBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
980e51e0e81SBarry Smith 
9813f9fe445SBarry Smith    Logically Collective on Mat
982fee21e36SBarry Smith 
983c7fcc2eaSBarry Smith     Input Parameters:
984c7fcc2eaSBarry Smith +   mat - the shell matrix
985c7fcc2eaSBarry Smith .   op - the name of the operation
986c7fcc2eaSBarry Smith -   f - the function that provides the operation.
987c7fcc2eaSBarry Smith 
98815091d37SBarry Smith    Level: advanced
98915091d37SBarry Smith 
990fae171e0SBarry Smith     Usage:
991e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
992f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
993c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
9940b627109SLois Curfman McInnes 
995a62d957aSLois Curfman McInnes     Notes:
996e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
9971c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
998a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
9991c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
1000a62d957aSLois Curfman McInnes 
1001a4d85fd7SBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
1002deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
1003deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
1004deebb3c3SLois Curfman McInnes     routines, e.g.,
1005a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1006a62d957aSLois Curfman McInnes 
10074aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
10084aa34b0aSBarry Smith     nonzero on failure.
10094aa34b0aSBarry Smith 
1010a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
1011a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
1012a62d957aSLois Curfman McInnes     set by MatCreateShell().
1013a62d957aSLois Curfman McInnes 
10142a7a6963SBarry Smith     Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
1015501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
1016501d9185SBarry Smith 
10170c0fd78eSBarry Smith     Use MatSetOperation() to set an operation for any matrix type
10180c0fd78eSBarry Smith 
1019a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
1020a62d957aSLois Curfman McInnes 
10210c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts()
1022e51e0e81SBarry Smith @*/
10237087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
1024e51e0e81SBarry Smith {
1025dfbe8321SBarry Smith   PetscErrorCode ierr;
1026ace3abfcSBarry Smith   PetscBool      flg;
10270c0fd78eSBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1028273d9f13SBarry Smith 
10293a40ed3dSBarry Smith   PetscFunctionBegin;
10300700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
10310c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
10320c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
10335edf6516SJed Brown   switch (op) {
10345edf6516SJed Brown   case MATOP_DESTROY:
103528f357e3SAlex Fikl     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
10365edf6516SJed Brown     break;
10376464bf51SAlex Fikl   case MATOP_DIAGONAL_SET:
10380c0fd78eSBarry Smith   case MATOP_DIAGONAL_SCALE:
10390c0fd78eSBarry Smith   case MATOP_SHIFT:
10400c0fd78eSBarry Smith   case MATOP_SCALE:
10419f137db4SBarry Smith   case MATOP_AXPY:
10420c0fd78eSBarry Smith     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
10430c0fd78eSBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
10446464bf51SAlex Fikl     break;
10450c0fd78eSBarry Smith   case MATOP_GET_DIAGONAL:
10460c0fd78eSBarry Smith     if (shell->managescalingshifts) shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
10470c0fd78eSBarry Smith     else mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
10485edf6516SJed Brown   case MATOP_VIEW:
104940a5a6c4SBarry Smith     if (!mat->ops->viewnative) {
105040a5a6c4SBarry Smith       mat->ops->viewnative = mat->ops->view;
105140a5a6c4SBarry Smith     }
10525edf6516SJed Brown     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
10535edf6516SJed Brown     break;
10545edf6516SJed Brown   case MATOP_MULT:
10559f137db4SBarry Smith     if (shell->managescalingshifts){
10569f137db4SBarry Smith       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10579f137db4SBarry Smith       mat->ops->mult   = MatMult_Shell;
10589f137db4SBarry Smith     } else mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10595edf6516SJed Brown     break;
10605edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
10619f137db4SBarry Smith     if (shell->managescalingshifts) {
10629f137db4SBarry Smith       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1063*5259c5a4SBarry Smith       mat->ops->multtranspose   = MatMultTranspose_Shell;
10649f137db4SBarry Smith     } else mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10655edf6516SJed Brown     break;
10665edf6516SJed Brown   default:
10675edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
1068a62d957aSLois Curfman McInnes   }
10693a40ed3dSBarry Smith   PetscFunctionReturn(0);
1070e51e0e81SBarry Smith }
1071f0479e8cSBarry Smith 
1072d4bb536fSBarry Smith /*@C
1073d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
1074d4bb536fSBarry Smith 
1075c7fcc2eaSBarry Smith     Not Collective
1076c7fcc2eaSBarry Smith 
1077d4bb536fSBarry Smith     Input Parameters:
1078c7fcc2eaSBarry Smith +   mat - the shell matrix
1079c7fcc2eaSBarry Smith -   op - the name of the operation
1080d4bb536fSBarry Smith 
1081d4bb536fSBarry Smith     Output Parameter:
1082d4bb536fSBarry Smith .   f - the function that provides the operation.
1083d4bb536fSBarry Smith 
108415091d37SBarry Smith     Level: advanced
108515091d37SBarry Smith 
1086d4bb536fSBarry Smith     Notes:
1087e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
1088d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
1089d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
1090d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
1091d4bb536fSBarry Smith 
1092d4bb536fSBarry Smith     All user-provided functions have the same calling
1093d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
1094d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
1095d4bb536fSBarry Smith     routines, e.g.,
1096d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1097d4bb536fSBarry Smith 
1098d4bb536fSBarry Smith     Within each user-defined routine, the user should call
1099d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
1100d4bb536fSBarry Smith     set by MatCreateShell().
1101d4bb536fSBarry Smith 
1102d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
1103d4bb536fSBarry Smith 
1104ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1105d4bb536fSBarry Smith @*/
11067087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
1107d4bb536fSBarry Smith {
1108dfbe8321SBarry Smith   PetscErrorCode ierr;
1109ace3abfcSBarry Smith   PetscBool      flg;
1110273d9f13SBarry Smith 
11113a40ed3dSBarry Smith   PetscFunctionBegin;
11120700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
111328f357e3SAlex Fikl   switch (op) {
111428f357e3SAlex Fikl   case MATOP_DESTROY:
1115251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1116273d9f13SBarry Smith     if (flg) {
1117d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
111828f357e3SAlex Fikl       *f = (void (*)(void))shell->ops->destroy;
1119c7fcc2eaSBarry Smith     } else {
1120c134de8dSSatish Balay       *f = (void (*)(void))mat->ops->destroy;
1121d4bb536fSBarry Smith     }
112228f357e3SAlex Fikl     break;
112328f357e3SAlex Fikl   case MATOP_DIAGONAL_SET:
112428f357e3SAlex Fikl     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
112528f357e3SAlex Fikl     if (flg) {
112628f357e3SAlex Fikl       Mat_Shell *shell = (Mat_Shell*)mat->data;
112728f357e3SAlex Fikl       *f = (void (*)(void))shell->ops->diagonalset;
1128c7fcc2eaSBarry Smith     } else {
112928f357e3SAlex Fikl       *f = (void (*)(void))mat->ops->diagonalset;
113028f357e3SAlex Fikl     }
113128f357e3SAlex Fikl     break;
113228f357e3SAlex Fikl   case MATOP_VIEW:
113328f357e3SAlex Fikl     *f = (void (*)(void))mat->ops->view;
113428f357e3SAlex Fikl     break;
113528f357e3SAlex Fikl   default:
1136c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
1137d4bb536fSBarry Smith   }
11383a40ed3dSBarry Smith   PetscFunctionReturn(0);
1139d4bb536fSBarry Smith }
1140d4bb536fSBarry Smith 
11410c0fd78eSBarry Smith /*@C
11420c0fd78eSBarry Smith     MatSetOperation - Allows user to set a matrix operation for any matrix type
11430c0fd78eSBarry Smith 
11440c0fd78eSBarry Smith    Logically Collective on Mat
11450c0fd78eSBarry Smith 
11460c0fd78eSBarry Smith     Input Parameters:
11470c0fd78eSBarry Smith +   mat - the shell matrix
11480c0fd78eSBarry Smith .   op - the name of the operation
11490c0fd78eSBarry Smith -   f - the function that provides the operation.
11500c0fd78eSBarry Smith 
11510c0fd78eSBarry Smith    Level: developer
11520c0fd78eSBarry Smith 
11530c0fd78eSBarry Smith     Usage:
11540c0fd78eSBarry Smith $      extern PetscErrorCode usermult(Mat,Vec,Vec);
11550c0fd78eSBarry Smith $      ierr = MatCreateXXX(comm,...&A);
11560c0fd78eSBarry Smith $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
11570c0fd78eSBarry Smith 
11580c0fd78eSBarry Smith     Notes:
11590c0fd78eSBarry Smith     See the file include/petscmat.h for a complete list of matrix
11600c0fd78eSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
11610c0fd78eSBarry Smith     <OPERATION> is the name (in all capital letters) of the
11620c0fd78eSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
11630c0fd78eSBarry Smith 
11640c0fd78eSBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
11650c0fd78eSBarry Smith     sequence as the usual matrix interface routines, since they
11660c0fd78eSBarry Smith     are intended to be accessed via the usual matrix interface
11670c0fd78eSBarry Smith     routines, e.g.,
11680c0fd78eSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
11690c0fd78eSBarry Smith 
11700c0fd78eSBarry Smith     In particular each function MUST return an error code of 0 on success and
11710c0fd78eSBarry Smith     nonzero on failure.
11720c0fd78eSBarry Smith 
11730c0fd78eSBarry Smith     This routine is distinct from MatShellSetOperation() in that it can be called on any matrix type.
11740c0fd78eSBarry Smith 
11750c0fd78eSBarry Smith .keywords: matrix, shell, set, operation
11760c0fd78eSBarry Smith 
11770c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
11780c0fd78eSBarry Smith @*/
11790c0fd78eSBarry Smith PetscErrorCode  MatSetOperation(Mat mat,MatOperation op,void (*f)(void))
11800c0fd78eSBarry Smith {
11810c0fd78eSBarry Smith   PetscFunctionBegin;
11820c0fd78eSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
11830c0fd78eSBarry Smith   (((void(**)(void))mat->ops)[op]) = f;
11840c0fd78eSBarry Smith   PetscFunctionReturn(0);
11850c0fd78eSBarry Smith }
1186