xref: /petsc/src/mat/impls/shell/shell.c (revision 81c0251941fee5a4fe0317982ad0ce67fac32af2)
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);}
376*81c02519SBarry Smith   if (shell->axpy) {
377*81c02519SBarry Smith     if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);}
378*81c02519SBarry Smith     ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr);
379*81c02519SBarry Smith     ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
380*81c02519SBarry Smith   }
38118be62a5SSatish Balay   PetscFunctionReturn(0);
38218be62a5SSatish Balay }
38318be62a5SSatish Balay 
384f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
385ef66eb69SBarry Smith {
386ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3878fe8eb27SJed Brown   PetscErrorCode ierr;
388b24ad042SBarry Smith 
389ef66eb69SBarry Smith   PetscFunctionBegin;
3900c0fd78eSBarry Smith   if (shell->left || shell->right) {
3918fe8eb27SJed Brown     if (!shell->dshift) {
3920c0fd78eSBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
3930c0fd78eSBarry Smith       ierr = VecSet(shell->dshift,a);CHKERRQ(ierr);
3940c0fd78eSBarry Smith     } else {
3950c0fd78eSBarry Smith       if (shell->left)  {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
3960c0fd78eSBarry Smith       if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
3979f137db4SBarry Smith       ierr = VecShift(shell->dshift,a);CHKERRQ(ierr);
3980c0fd78eSBarry Smith     }
3998fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
4008fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
4018fe8eb27SJed Brown   } else shell->vshift += a;
402ef66eb69SBarry Smith   PetscFunctionReturn(0);
403ef66eb69SBarry Smith }
404ef66eb69SBarry Smith 
4056464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
4066464bf51SAlex Fikl {
4076464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
4086464bf51SAlex Fikl   PetscErrorCode ierr;
4096464bf51SAlex Fikl 
4106464bf51SAlex Fikl   PetscFunctionBegin;
4110c0fd78eSBarry Smith   if (ins == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported with INSERT_VALUES");
4120c0fd78eSBarry Smith   if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);}
4130c0fd78eSBarry Smith   if (shell->left || shell->right) {
4149f137db4SBarry Smith     if (!shell->right_work) ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);
4159f137db4SBarry Smith     if (shell->left && shell->right)  {
4169f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
4179f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
4189f137db4SBarry Smith     } else if (shell->left) {
4199f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
4209f137db4SBarry Smith     } else {
4219f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr);
4229f137db4SBarry Smith     }
4239f137db4SBarry Smith     if (!shell->dshift) {
4249f137db4SBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
4259f137db4SBarry Smith       ierr = VecCopy(shell->dshift,shell->right_work);CHKERRQ(ierr);
4269f137db4SBarry Smith     } else {
4279f137db4SBarry Smith       ierr = VecAXPY(shell->dshift,1.0,shell->right_work);CHKERRQ(ierr);
4289f137db4SBarry Smith     }
4290c0fd78eSBarry Smith   } else {
4300c0fd78eSBarry Smith     ierr = VecAXPY(shell->dshift,1.0,D);CHKERRQ(ierr);
4316464bf51SAlex Fikl   }
4326464bf51SAlex Fikl   PetscFunctionReturn(0);
4336464bf51SAlex Fikl }
4346464bf51SAlex Fikl 
435f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
436ef66eb69SBarry Smith {
437ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4388fe8eb27SJed Brown   PetscErrorCode ierr;
439b24ad042SBarry Smith 
440ef66eb69SBarry Smith   PetscFunctionBegin;
441f4df32b1SMatthew Knepley   shell->vscale *= a;
4420c0fd78eSBarry Smith   shell->vshift *= a;
4432205254eSKarl Rupp   if (shell->dshift) {
4442205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
4450c0fd78eSBarry Smith   }
446*81c02519SBarry Smith   shell->axpy_vscale *= a;
4478fe8eb27SJed Brown   PetscFunctionReturn(0);
44818be62a5SSatish Balay }
4498fe8eb27SJed Brown 
4508fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
4518fe8eb27SJed Brown {
4528fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4538fe8eb27SJed Brown   PetscErrorCode ierr;
4548fe8eb27SJed Brown 
4558fe8eb27SJed Brown   PetscFunctionBegin;
456*81c02519SBarry Smith   if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation");
4578fe8eb27SJed Brown   if (left) {
4580c0fd78eSBarry Smith     if (!shell->left) {
4590c0fd78eSBarry Smith       ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr);
4608fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
4610c0fd78eSBarry Smith     } else {
4620c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
46318be62a5SSatish Balay     }
464ef66eb69SBarry Smith   }
4658fe8eb27SJed Brown   if (right) {
4660c0fd78eSBarry Smith     if (!shell->right) {
4670c0fd78eSBarry Smith       ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr);
4688fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
4690c0fd78eSBarry Smith     } else {
4700c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
4718fe8eb27SJed Brown     }
4728fe8eb27SJed Brown   }
473ef66eb69SBarry Smith   PetscFunctionReturn(0);
474ef66eb69SBarry Smith }
475ef66eb69SBarry Smith 
476dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
477ef66eb69SBarry Smith {
478ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4790c0fd78eSBarry Smith   PetscErrorCode ierr;
480ef66eb69SBarry Smith 
481ef66eb69SBarry Smith   PetscFunctionBegin;
4828fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
483ef66eb69SBarry Smith     shell->vshift = 0.0;
484ef66eb69SBarry Smith     shell->vscale = 1.0;
4850c0fd78eSBarry Smith     ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
4860c0fd78eSBarry Smith     ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
4870c0fd78eSBarry Smith     ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
488*81c02519SBarry Smith     ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
489ef66eb69SBarry Smith   }
490ef66eb69SBarry Smith   PetscFunctionReturn(0);
491ef66eb69SBarry Smith }
492ef66eb69SBarry Smith 
493cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
494b951964fSBarry Smith 
4953b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
4963b49f96aSBarry Smith {
4973b49f96aSBarry Smith   PetscFunctionBegin;
4983b49f96aSBarry Smith   *missing = PETSC_FALSE;
4993b49f96aSBarry Smith   PetscFunctionReturn(0);
5003b49f96aSBarry Smith }
5013b49f96aSBarry Smith 
5029f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
5039f137db4SBarry Smith {
5049f137db4SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
5059f137db4SBarry Smith   PetscErrorCode ierr;
5069f137db4SBarry Smith 
5079f137db4SBarry Smith   PetscFunctionBegin;
5089f137db4SBarry Smith   ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr);
5099f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
5109f137db4SBarry Smith   shell->axpy        = X;
5119f137db4SBarry Smith   shell->axpy_vscale = a;
5129f137db4SBarry Smith   PetscFunctionReturn(0);
5139f137db4SBarry Smith }
5149f137db4SBarry Smith 
51509dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
51620563c6bSBarry Smith                                        0,
51720563c6bSBarry Smith                                        0,
5189f137db4SBarry Smith                                        0,
5190c0fd78eSBarry Smith                                 /* 4*/ MatMultAdd_Shell,
5209f137db4SBarry Smith                                        0,
5210c0fd78eSBarry Smith                                        MatMultTransposeAdd_Shell,
522b951964fSBarry Smith                                        0,
523b951964fSBarry Smith                                        0,
524b951964fSBarry Smith                                        0,
52597304618SKris Buschelman                                 /*10*/ 0,
526b951964fSBarry Smith                                        0,
527b951964fSBarry Smith                                        0,
528b951964fSBarry Smith                                        0,
529b951964fSBarry Smith                                        0,
53097304618SKris Buschelman                                 /*15*/ 0,
531b951964fSBarry Smith                                        0,
5320c0fd78eSBarry Smith                                        MatGetDiagonal_Shell,
5338fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
534b951964fSBarry Smith                                        0,
53597304618SKris Buschelman                                 /*20*/ 0,
536ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
537b951964fSBarry Smith                                        0,
538b951964fSBarry Smith                                        0,
539d519adbfSMatthew Knepley                                 /*24*/ 0,
540b951964fSBarry Smith                                        0,
541b951964fSBarry Smith                                        0,
542b951964fSBarry Smith                                        0,
543b951964fSBarry Smith                                        0,
544d519adbfSMatthew Knepley                                 /*29*/ 0,
545b951964fSBarry Smith                                        0,
546273d9f13SBarry Smith                                        0,
547b951964fSBarry Smith                                        0,
548b951964fSBarry Smith                                        0,
549cb8c8a77SBarry Smith                                 /*34*/ MatDuplicate_Shell,
550b951964fSBarry Smith                                        0,
551b951964fSBarry Smith                                        0,
55209dc0095SBarry Smith                                        0,
55309dc0095SBarry Smith                                        0,
5549f137db4SBarry Smith                                 /*39*/ MatAXPY_Shell,
55509dc0095SBarry Smith                                        0,
55609dc0095SBarry Smith                                        0,
55709dc0095SBarry Smith                                        0,
558cb8c8a77SBarry Smith                                        MatCopy_Shell,
559d519adbfSMatthew Knepley                                 /*44*/ 0,
560ef66eb69SBarry Smith                                        MatScale_Shell,
561ef66eb69SBarry Smith                                        MatShift_Shell,
5626464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
56309dc0095SBarry Smith                                        0,
564f73d5cc4SBarry Smith                                 /*49*/ 0,
56509dc0095SBarry Smith                                        0,
56609dc0095SBarry Smith                                        0,
56709dc0095SBarry Smith                                        0,
56809dc0095SBarry Smith                                        0,
569d519adbfSMatthew Knepley                                 /*54*/ 0,
57009dc0095SBarry Smith                                        0,
57109dc0095SBarry Smith                                        0,
57209dc0095SBarry Smith                                        0,
57309dc0095SBarry Smith                                        0,
574d519adbfSMatthew Knepley                                 /*59*/ 0,
575b9b97703SBarry Smith                                        MatDestroy_Shell,
57609dc0095SBarry Smith                                        0,
577357abbc8SBarry Smith                                        0,
578273d9f13SBarry Smith                                        0,
579d519adbfSMatthew Knepley                                 /*64*/ 0,
580273d9f13SBarry Smith                                        0,
581273d9f13SBarry Smith                                        0,
582273d9f13SBarry Smith                                        0,
583273d9f13SBarry Smith                                        0,
584d519adbfSMatthew Knepley                                 /*69*/ 0,
585273d9f13SBarry Smith                                        0,
586c87e5d42SMatthew Knepley                                        MatConvert_Shell,
587273d9f13SBarry Smith                                        0,
58897304618SKris Buschelman                                        0,
589d519adbfSMatthew Knepley                                 /*74*/ 0,
59097304618SKris Buschelman                                        0,
59197304618SKris Buschelman                                        0,
59297304618SKris Buschelman                                        0,
59397304618SKris Buschelman                                        0,
594d519adbfSMatthew Knepley                                 /*79*/ 0,
59597304618SKris Buschelman                                        0,
59697304618SKris Buschelman                                        0,
59797304618SKris Buschelman                                        0,
598865e5f61SKris Buschelman                                        0,
599d519adbfSMatthew Knepley                                 /*84*/ 0,
600865e5f61SKris Buschelman                                        0,
601865e5f61SKris Buschelman                                        0,
602865e5f61SKris Buschelman                                        0,
603865e5f61SKris Buschelman                                        0,
604d519adbfSMatthew Knepley                                 /*89*/ 0,
605865e5f61SKris Buschelman                                        0,
606865e5f61SKris Buschelman                                        0,
607865e5f61SKris Buschelman                                        0,
608865e5f61SKris Buschelman                                        0,
609d519adbfSMatthew Knepley                                 /*94*/ 0,
610865e5f61SKris Buschelman                                        0,
611865e5f61SKris Buschelman                                        0,
6123964eb88SJed Brown                                        0,
6133964eb88SJed Brown                                        0,
6143964eb88SJed Brown                                 /*99*/ 0,
6153964eb88SJed Brown                                        0,
6163964eb88SJed Brown                                        0,
6173964eb88SJed Brown                                        0,
6183964eb88SJed Brown                                        0,
6193964eb88SJed Brown                                /*104*/ 0,
6203964eb88SJed Brown                                        0,
6213964eb88SJed Brown                                        0,
6223964eb88SJed Brown                                        0,
6233964eb88SJed Brown                                        0,
6243964eb88SJed Brown                                /*109*/ 0,
6253964eb88SJed Brown                                        0,
6263964eb88SJed Brown                                        0,
6273964eb88SJed Brown                                        0,
6283b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
6293964eb88SJed Brown                                /*114*/ 0,
6303964eb88SJed Brown                                        0,
6313964eb88SJed Brown                                        0,
6323964eb88SJed Brown                                        0,
6333964eb88SJed Brown                                        0,
6343964eb88SJed Brown                                /*119*/ 0,
6353964eb88SJed Brown                                        0,
6363964eb88SJed Brown                                        0,
6373964eb88SJed Brown                                        0,
6383964eb88SJed Brown                                        0,
6393964eb88SJed Brown                                /*124*/ 0,
6403964eb88SJed Brown                                        0,
6413964eb88SJed Brown                                        0,
6423964eb88SJed Brown                                        0,
6433964eb88SJed Brown                                        0,
6443964eb88SJed Brown                                /*129*/ 0,
6453964eb88SJed Brown                                        0,
6463964eb88SJed Brown                                        0,
6473964eb88SJed Brown                                        0,
6483964eb88SJed Brown                                        0,
6493964eb88SJed Brown                                /*134*/ 0,
6503964eb88SJed Brown                                        0,
6513964eb88SJed Brown                                        0,
6523964eb88SJed Brown                                        0,
6533964eb88SJed Brown                                        0,
6543964eb88SJed Brown                                /*139*/ 0,
655f9426fe0SMark Adams                                        0,
6563964eb88SJed Brown                                        0
6573964eb88SJed Brown };
658273d9f13SBarry Smith 
6590bad9183SKris Buschelman /*MC
660fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
6610bad9183SKris Buschelman 
6620bad9183SKris Buschelman   Level: advanced
6630bad9183SKris Buschelman 
6640c0fd78eSBarry Smith .seealso: MatCreateShell()
6650bad9183SKris Buschelman M*/
6660bad9183SKris Buschelman 
6678cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
668273d9f13SBarry Smith {
669273d9f13SBarry Smith   Mat_Shell      *b;
670dfbe8321SBarry Smith   PetscErrorCode ierr;
671273d9f13SBarry Smith 
672273d9f13SBarry Smith   PetscFunctionBegin;
673273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
674273d9f13SBarry Smith 
675b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
676273d9f13SBarry Smith   A->data = (void*)b;
677273d9f13SBarry Smith 
67826283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
67926283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
680273d9f13SBarry Smith 
681273d9f13SBarry Smith   b->ctx                 = 0;
682ef66eb69SBarry Smith   b->vshift              = 0.0;
683ef66eb69SBarry Smith   b->vscale              = 1.0;
6840c0fd78eSBarry Smith   b->managescalingshifts = PETSC_TRUE;
685273d9f13SBarry Smith   A->assembled           = PETSC_TRUE;
686210f0121SBarry Smith   A->preallocated        = PETSC_FALSE;
6872205254eSKarl Rupp 
68817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
689273d9f13SBarry Smith   PetscFunctionReturn(0);
690273d9f13SBarry Smith }
691e51e0e81SBarry Smith 
6924b828684SBarry Smith /*@C
693052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
694ff756334SLois Curfman McInnes    private data storage format.
695e51e0e81SBarry Smith 
696c7fcc2eaSBarry Smith   Collective on MPI_Comm
697c7fcc2eaSBarry Smith 
698e51e0e81SBarry Smith    Input Parameters:
699c7fcc2eaSBarry Smith +  comm - MPI communicator
700c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
701c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
702c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
703c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
704c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
705e51e0e81SBarry Smith 
706ff756334SLois Curfman McInnes    Output Parameter:
70744cd7ae7SLois Curfman McInnes .  A - the matrix
708e51e0e81SBarry Smith 
709ff2fd236SBarry Smith    Level: advanced
710ff2fd236SBarry Smith 
711f39d1f56SLois Curfman McInnes   Usage:
7127b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
713f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
714c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
715f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
716f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
717f39d1f56SLois Curfman McInnes 
718ff756334SLois Curfman McInnes    Notes:
719ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
720ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
721ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
722e51e0e81SBarry Smith 
723daf670e6SBarry Smith    Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this
724daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
725daf670e6SBarry Smith     in as the ctx argument.
726f38a8d56SBarry Smith 
727f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
728f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
729645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
730645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
731645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
732645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
733645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
734645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
735645985a0SLois Curfman McInnes    For example,
736f39d1f56SLois Curfman McInnes 
7370c0fd78eSBarry Smith    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), and MatScale() internally so these
7380c0fd78eSBarry Smith    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.
7390c0fd78eSBarry Smith 
740f39d1f56SLois Curfman McInnes $
741f39d1f56SLois Curfman McInnes $     Vec x, y
7427b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
743645985a0SLois Curfman McInnes $     Mat A
744f39d1f56SLois Curfman McInnes $
745c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
746c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
747f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
748c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
749c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
750c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
751645985a0SLois Curfman McInnes $     MatMult(A,x,y);
752645985a0SLois Curfman McInnes $     MatDestroy(A);
753f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
754645985a0SLois Curfman McInnes $
755e51e0e81SBarry Smith 
7560c0fd78eSBarry Smith     For rectangular matrices do all the scalings and shifts make sense?
7570c0fd78eSBarry Smith 
7580c0fd78eSBarry Smith     Developers Notes: Regarding shifting and scaling. The general form is
7590c0fd78eSBarry Smith 
7600c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
7610c0fd78eSBarry Smith 
7620c0fd78eSBarry Smith       The order you apply the operations is important. For example if you have a dshift then
7630c0fd78eSBarry Smith       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
7640c0fd78eSBarry Smith       you get s*vscale*A + diag(shift)
7650c0fd78eSBarry Smith 
7660c0fd78eSBarry Smith           A is the user provided function.
7670c0fd78eSBarry Smith 
7680b627109SLois Curfman McInnes .keywords: matrix, shell, create
7690b627109SLois Curfman McInnes 
7700c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
771e51e0e81SBarry Smith @*/
7727087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
773e51e0e81SBarry Smith {
774dfbe8321SBarry Smith   PetscErrorCode ierr;
775ed3cc1f0SBarry Smith 
7763a40ed3dSBarry Smith   PetscFunctionBegin;
777f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
778f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
779273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
780273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
7810fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
782273d9f13SBarry Smith   PetscFunctionReturn(0);
783c7fcc2eaSBarry Smith }
784c7fcc2eaSBarry Smith 
785c6866cfdSSatish Balay /*@
786273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
787c7fcc2eaSBarry Smith 
7883f9fe445SBarry Smith    Logically Collective on Mat
789c7fcc2eaSBarry Smith 
790273d9f13SBarry Smith     Input Parameters:
791273d9f13SBarry Smith +   mat - the shell matrix
792273d9f13SBarry Smith -   ctx - the context
793273d9f13SBarry Smith 
794273d9f13SBarry Smith    Level: advanced
795273d9f13SBarry Smith 
796daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
797daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
798273d9f13SBarry Smith 
799273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
8000bc0a719SBarry Smith @*/
8017087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
802273d9f13SBarry Smith {
803273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
804dfbe8321SBarry Smith   PetscErrorCode ierr;
805ace3abfcSBarry Smith   PetscBool      flg;
806273d9f13SBarry Smith 
807273d9f13SBarry Smith   PetscFunctionBegin;
8080700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
809251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
810273d9f13SBarry Smith   if (flg) {
811273d9f13SBarry Smith     shell->ctx = ctx;
812ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
8133a40ed3dSBarry Smith   PetscFunctionReturn(0);
814e51e0e81SBarry Smith }
815e51e0e81SBarry Smith 
8160c0fd78eSBarry Smith /*@
8170c0fd78eSBarry Smith     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
8180c0fd78eSBarry Smith           after MatCreateShell()
8190c0fd78eSBarry Smith 
8200c0fd78eSBarry Smith    Logically Collective on Mat
8210c0fd78eSBarry Smith 
8220c0fd78eSBarry Smith     Input Parameter:
8230c0fd78eSBarry Smith .   mat - the shell matrix
8240c0fd78eSBarry Smith 
8250c0fd78eSBarry Smith   Level: advanced
8260c0fd78eSBarry Smith 
8270c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
8280c0fd78eSBarry Smith @*/
8290c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A)
8300c0fd78eSBarry Smith {
8310c0fd78eSBarry Smith   PetscErrorCode ierr;
8320c0fd78eSBarry Smith   Mat_Shell      *shell = (Mat_Shell*)A->data;
8330c0fd78eSBarry Smith   PetscBool      flg;
8340c0fd78eSBarry Smith 
8350c0fd78eSBarry Smith   PetscFunctionBegin;
8360c0fd78eSBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
8370c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr);
8380c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
8390c0fd78eSBarry Smith   shell->managescalingshifts = PETSC_FALSE;
8400c0fd78eSBarry Smith   A->ops->diagonalscale = 0;
8410c0fd78eSBarry Smith   A->ops->scale         = 0;
8420c0fd78eSBarry Smith   A->ops->shift         = 0;
8430c0fd78eSBarry Smith   A->ops->diagonalset   = 0;
8440c0fd78eSBarry Smith   PetscFunctionReturn(0);
8450c0fd78eSBarry Smith }
8460c0fd78eSBarry Smith 
847c16cb8f2SBarry Smith /*@C
848f3b1f45cSBarry Smith     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
849f3b1f45cSBarry Smith 
850f3b1f45cSBarry Smith    Logically Collective on Mat
851f3b1f45cSBarry Smith 
852f3b1f45cSBarry Smith     Input Parameters:
853f3b1f45cSBarry Smith +   mat - the shell matrix
854f3b1f45cSBarry Smith .   f - the function
855f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
856f3b1f45cSBarry Smith -   ctx - an optional context for the function
857f3b1f45cSBarry Smith 
858f3b1f45cSBarry Smith    Output Parameter:
859f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
860f3b1f45cSBarry Smith 
861f3b1f45cSBarry Smith    Options Database:
862f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
863f3b1f45cSBarry Smith 
864f3b1f45cSBarry Smith    Level: advanced
865f3b1f45cSBarry Smith 
866f3b1f45cSBarry Smith    Fortran Notes: Not supported from Fortran
867f3b1f45cSBarry Smith 
868f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
869f3b1f45cSBarry Smith @*/
870f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
871f3b1f45cSBarry Smith {
872f3b1f45cSBarry Smith   PetscErrorCode ierr;
873f3b1f45cSBarry Smith   PetscInt       m,n;
874f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
875f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
87674e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
877f3b1f45cSBarry Smith 
878f3b1f45cSBarry Smith   PetscFunctionBegin;
879f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
880f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr);
881f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
882f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr);
883f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
884f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
885f3b1f45cSBarry Smith 
886f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
887f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mat,&Dmat);CHKERRQ(ierr);
888f3b1f45cSBarry Smith 
889f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
890f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
891f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
892f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
893f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
894f3b1f45cSBarry Smith     flag = PETSC_FALSE;
895f3b1f45cSBarry Smith     if (v) {
896f3b1f45cSBarry 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));
897f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
898f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
899f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
900f3b1f45cSBarry Smith     }
901f3b1f45cSBarry Smith   } else if (v) {
902f3b1f45cSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");
903f3b1f45cSBarry Smith   }
904f3b1f45cSBarry Smith   if (flg) *flg = flag;
905f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
906f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
907f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
908f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
909f3b1f45cSBarry Smith   PetscFunctionReturn(0);
910f3b1f45cSBarry Smith }
911f3b1f45cSBarry Smith 
912f3b1f45cSBarry Smith /*@C
913f3b1f45cSBarry Smith     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.
914f3b1f45cSBarry Smith 
915f3b1f45cSBarry Smith    Logically Collective on Mat
916f3b1f45cSBarry Smith 
917f3b1f45cSBarry Smith     Input Parameters:
918f3b1f45cSBarry Smith +   mat - the shell matrix
919f3b1f45cSBarry Smith .   f - the function
920f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
921f3b1f45cSBarry Smith -   ctx - an optional context for the function
922f3b1f45cSBarry Smith 
923f3b1f45cSBarry Smith    Output Parameter:
924f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
925f3b1f45cSBarry Smith 
926f3b1f45cSBarry Smith    Options Database:
927f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
928f3b1f45cSBarry Smith 
929f3b1f45cSBarry Smith    Level: advanced
930f3b1f45cSBarry Smith 
931f3b1f45cSBarry Smith    Fortran Notes: Not supported from Fortran
932f3b1f45cSBarry Smith 
933f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
934f3b1f45cSBarry Smith @*/
935f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
936f3b1f45cSBarry Smith {
937f3b1f45cSBarry Smith   PetscErrorCode ierr;
938f3b1f45cSBarry Smith   Vec            x,y,z;
939f3b1f45cSBarry Smith   PetscInt       m,n,M,N;
940f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
941f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
94274e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
943f3b1f45cSBarry Smith 
944f3b1f45cSBarry Smith   PetscFunctionBegin;
945f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
946f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
947f3b1f45cSBarry Smith   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
948f3b1f45cSBarry Smith   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
949f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
950f3b1f45cSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
951f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
952f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
953f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
954f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
955f3b1f45cSBarry Smith   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
956f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr);
957f3b1f45cSBarry Smith 
958f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
959f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
960f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
961f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
962f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
963f3b1f45cSBarry Smith     flag = PETSC_FALSE;
964f3b1f45cSBarry Smith     if (v) {
965f3b1f45cSBarry 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));
966f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
967f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
968f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
969f3b1f45cSBarry Smith     }
970f3b1f45cSBarry Smith   } else if (v) {
971f3b1f45cSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");
972f3b1f45cSBarry Smith   }
973f3b1f45cSBarry Smith   if (flg) *flg = flag;
974f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
975f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
976f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
977f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
978f3b1f45cSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
979f3b1f45cSBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
980f3b1f45cSBarry Smith   ierr = VecDestroy(&z);CHKERRQ(ierr);
981f3b1f45cSBarry Smith   PetscFunctionReturn(0);
982f3b1f45cSBarry Smith }
983f3b1f45cSBarry Smith 
984f3b1f45cSBarry Smith /*@C
9853a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
9863a3eedf2SBarry Smith                            a shell matrix.
9870c0fd78eSBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
988e51e0e81SBarry Smith 
9893f9fe445SBarry Smith    Logically Collective on Mat
990fee21e36SBarry Smith 
991c7fcc2eaSBarry Smith     Input Parameters:
992c7fcc2eaSBarry Smith +   mat - the shell matrix
993c7fcc2eaSBarry Smith .   op - the name of the operation
994c7fcc2eaSBarry Smith -   f - the function that provides the operation.
995c7fcc2eaSBarry Smith 
99615091d37SBarry Smith    Level: advanced
99715091d37SBarry Smith 
998fae171e0SBarry Smith     Usage:
999e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
1000f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
1001c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
10020b627109SLois Curfman McInnes 
1003a62d957aSLois Curfman McInnes     Notes:
1004e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
10051c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
1006a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
10071c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
1008a62d957aSLois Curfman McInnes 
1009a4d85fd7SBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
1010deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
1011deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
1012deebb3c3SLois Curfman McInnes     routines, e.g.,
1013a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1014a62d957aSLois Curfman McInnes 
10154aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
10164aa34b0aSBarry Smith     nonzero on failure.
10174aa34b0aSBarry Smith 
1018a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
1019a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
1020a62d957aSLois Curfman McInnes     set by MatCreateShell().
1021a62d957aSLois Curfman McInnes 
10222a7a6963SBarry Smith     Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
1023501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
1024501d9185SBarry Smith 
10250c0fd78eSBarry Smith     Use MatSetOperation() to set an operation for any matrix type
10260c0fd78eSBarry Smith 
1027a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
1028a62d957aSLois Curfman McInnes 
10290c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts()
1030e51e0e81SBarry Smith @*/
10317087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
1032e51e0e81SBarry Smith {
1033dfbe8321SBarry Smith   PetscErrorCode ierr;
1034ace3abfcSBarry Smith   PetscBool      flg;
10350c0fd78eSBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1036273d9f13SBarry Smith 
10373a40ed3dSBarry Smith   PetscFunctionBegin;
10380700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
10390c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
10400c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
10415edf6516SJed Brown   switch (op) {
10425edf6516SJed Brown   case MATOP_DESTROY:
104328f357e3SAlex Fikl     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
10445edf6516SJed Brown     break;
10456464bf51SAlex Fikl   case MATOP_DIAGONAL_SET:
10460c0fd78eSBarry Smith   case MATOP_DIAGONAL_SCALE:
10470c0fd78eSBarry Smith   case MATOP_SHIFT:
10480c0fd78eSBarry Smith   case MATOP_SCALE:
10499f137db4SBarry Smith   case MATOP_AXPY:
10500c0fd78eSBarry Smith     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
10510c0fd78eSBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
10526464bf51SAlex Fikl     break;
10530c0fd78eSBarry Smith   case MATOP_GET_DIAGONAL:
10540c0fd78eSBarry Smith     if (shell->managescalingshifts) shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
10550c0fd78eSBarry Smith     else mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
10565edf6516SJed Brown   case MATOP_VIEW:
105740a5a6c4SBarry Smith     if (!mat->ops->viewnative) {
105840a5a6c4SBarry Smith       mat->ops->viewnative = mat->ops->view;
105940a5a6c4SBarry Smith     }
10605edf6516SJed Brown     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
10615edf6516SJed Brown     break;
10625edf6516SJed Brown   case MATOP_MULT:
10639f137db4SBarry Smith     if (shell->managescalingshifts){
10649f137db4SBarry Smith       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10659f137db4SBarry Smith       mat->ops->mult   = MatMult_Shell;
10669f137db4SBarry Smith     } else mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10675edf6516SJed Brown     break;
10685edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
10699f137db4SBarry Smith     if (shell->managescalingshifts) {
10709f137db4SBarry Smith       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10715259c5a4SBarry Smith       mat->ops->multtranspose   = MatMultTranspose_Shell;
10729f137db4SBarry Smith     } else mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10735edf6516SJed Brown     break;
10745edf6516SJed Brown   default:
10755edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
1076a62d957aSLois Curfman McInnes   }
10773a40ed3dSBarry Smith   PetscFunctionReturn(0);
1078e51e0e81SBarry Smith }
1079f0479e8cSBarry Smith 
1080d4bb536fSBarry Smith /*@C
1081d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
1082d4bb536fSBarry Smith 
1083c7fcc2eaSBarry Smith     Not Collective
1084c7fcc2eaSBarry Smith 
1085d4bb536fSBarry Smith     Input Parameters:
1086c7fcc2eaSBarry Smith +   mat - the shell matrix
1087c7fcc2eaSBarry Smith -   op - the name of the operation
1088d4bb536fSBarry Smith 
1089d4bb536fSBarry Smith     Output Parameter:
1090d4bb536fSBarry Smith .   f - the function that provides the operation.
1091d4bb536fSBarry Smith 
109215091d37SBarry Smith     Level: advanced
109315091d37SBarry Smith 
1094d4bb536fSBarry Smith     Notes:
1095e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
1096d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
1097d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
1098d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
1099d4bb536fSBarry Smith 
1100d4bb536fSBarry Smith     All user-provided functions have the same calling
1101d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
1102d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
1103d4bb536fSBarry Smith     routines, e.g.,
1104d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1105d4bb536fSBarry Smith 
1106d4bb536fSBarry Smith     Within each user-defined routine, the user should call
1107d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
1108d4bb536fSBarry Smith     set by MatCreateShell().
1109d4bb536fSBarry Smith 
1110d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
1111d4bb536fSBarry Smith 
1112ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1113d4bb536fSBarry Smith @*/
11147087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
1115d4bb536fSBarry Smith {
1116dfbe8321SBarry Smith   PetscErrorCode ierr;
1117ace3abfcSBarry Smith   PetscBool      flg;
1118273d9f13SBarry Smith 
11193a40ed3dSBarry Smith   PetscFunctionBegin;
11200700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
112128f357e3SAlex Fikl   switch (op) {
112228f357e3SAlex Fikl   case MATOP_DESTROY:
1123251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1124273d9f13SBarry Smith     if (flg) {
1125d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
112628f357e3SAlex Fikl       *f = (void (*)(void))shell->ops->destroy;
1127c7fcc2eaSBarry Smith     } else {
1128c134de8dSSatish Balay       *f = (void (*)(void))mat->ops->destroy;
1129d4bb536fSBarry Smith     }
113028f357e3SAlex Fikl     break;
113128f357e3SAlex Fikl   case MATOP_DIAGONAL_SET:
113228f357e3SAlex Fikl     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
113328f357e3SAlex Fikl     if (flg) {
113428f357e3SAlex Fikl       Mat_Shell *shell = (Mat_Shell*)mat->data;
113528f357e3SAlex Fikl       *f = (void (*)(void))shell->ops->diagonalset;
1136c7fcc2eaSBarry Smith     } else {
113728f357e3SAlex Fikl       *f = (void (*)(void))mat->ops->diagonalset;
113828f357e3SAlex Fikl     }
113928f357e3SAlex Fikl     break;
114028f357e3SAlex Fikl   case MATOP_VIEW:
114128f357e3SAlex Fikl     *f = (void (*)(void))mat->ops->view;
114228f357e3SAlex Fikl     break;
114328f357e3SAlex Fikl   default:
1144c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
1145d4bb536fSBarry Smith   }
11463a40ed3dSBarry Smith   PetscFunctionReturn(0);
1147d4bb536fSBarry Smith }
1148d4bb536fSBarry Smith 
11490c0fd78eSBarry Smith /*@C
11500c0fd78eSBarry Smith     MatSetOperation - Allows user to set a matrix operation for any matrix type
11510c0fd78eSBarry Smith 
11520c0fd78eSBarry Smith    Logically Collective on Mat
11530c0fd78eSBarry Smith 
11540c0fd78eSBarry Smith     Input Parameters:
11550c0fd78eSBarry Smith +   mat - the shell matrix
11560c0fd78eSBarry Smith .   op - the name of the operation
11570c0fd78eSBarry Smith -   f - the function that provides the operation.
11580c0fd78eSBarry Smith 
11590c0fd78eSBarry Smith    Level: developer
11600c0fd78eSBarry Smith 
11610c0fd78eSBarry Smith     Usage:
11620c0fd78eSBarry Smith $      extern PetscErrorCode usermult(Mat,Vec,Vec);
11630c0fd78eSBarry Smith $      ierr = MatCreateXXX(comm,...&A);
11640c0fd78eSBarry Smith $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
11650c0fd78eSBarry Smith 
11660c0fd78eSBarry Smith     Notes:
11670c0fd78eSBarry Smith     See the file include/petscmat.h for a complete list of matrix
11680c0fd78eSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
11690c0fd78eSBarry Smith     <OPERATION> is the name (in all capital letters) of the
11700c0fd78eSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
11710c0fd78eSBarry Smith 
11720c0fd78eSBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
11730c0fd78eSBarry Smith     sequence as the usual matrix interface routines, since they
11740c0fd78eSBarry Smith     are intended to be accessed via the usual matrix interface
11750c0fd78eSBarry Smith     routines, e.g.,
11760c0fd78eSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
11770c0fd78eSBarry Smith 
11780c0fd78eSBarry Smith     In particular each function MUST return an error code of 0 on success and
11790c0fd78eSBarry Smith     nonzero on failure.
11800c0fd78eSBarry Smith 
11810c0fd78eSBarry Smith     This routine is distinct from MatShellSetOperation() in that it can be called on any matrix type.
11820c0fd78eSBarry Smith 
11830c0fd78eSBarry Smith .keywords: matrix, shell, set, operation
11840c0fd78eSBarry Smith 
11850c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
11860c0fd78eSBarry Smith @*/
11870c0fd78eSBarry Smith PetscErrorCode  MatSetOperation(Mat mat,MatOperation op,void (*f)(void))
11880c0fd78eSBarry Smith {
11890c0fd78eSBarry Smith   PetscFunctionBegin;
11900c0fd78eSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
11910c0fd78eSBarry Smith   (((void(**)(void))mat->ops)[op]) = f;
11920c0fd78eSBarry Smith   PetscFunctionReturn(0);
11930c0fd78eSBarry Smith }
1194