xref: /petsc/src/mat/impls/shell/shell.c (revision b253ae0bb0f47342fc3f41d41f137728d93ccdb1)
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 {
11976e8c5aSLisandro Dalcin   /*  3 */ PetscErrorCode (*mult)(Mat,Vec,Vec);
12976e8c5aSLisandro Dalcin   /*  5 */ PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
13976e8c5aSLisandro Dalcin   /* 17 */ PetscErrorCode (*getdiagonal)(Mat,Vec);
14976e8c5aSLisandro Dalcin   /* 43 */ PetscErrorCode (*copy)(Mat,Mat,MatStructure);
15976e8c5aSLisandro Dalcin   /* 60 */ PetscErrorCode (*destroy)(Mat);
1628f357e3SAlex Fikl };
1728f357e3SAlex Fikl 
1828f357e3SAlex Fikl typedef struct {
1928f357e3SAlex Fikl   struct _MatShellOps ops[1];
202205254eSKarl Rupp 
21ef66eb69SBarry Smith   PetscScalar vscale,vshift;
228fe8eb27SJed Brown   Vec         dshift;
238fe8eb27SJed Brown   Vec         left,right;
248fe8eb27SJed Brown   Vec         left_work,right_work;
255edf6516SJed Brown   Vec         left_add_work,right_add_work;
269f137db4SBarry Smith   Mat         axpy;
279f137db4SBarry Smith   PetscScalar axpy_vscale;
280c0fd78eSBarry Smith   PetscBool   managescalingshifts;                   /* The user will manage the scaling and shifts for the MATSHELL, not the default */
2920563c6bSBarry Smith   void        *ctx;
3088cf3e7dSBarry Smith } Mat_Shell;
310c0fd78eSBarry Smith 
328fe8eb27SJed Brown /*
330c0fd78eSBarry Smith       xx = diag(left)*x
348fe8eb27SJed Brown */
358fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
368fe8eb27SJed Brown {
378fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
388fe8eb27SJed Brown   PetscErrorCode ierr;
398fe8eb27SJed Brown 
408fe8eb27SJed Brown   PetscFunctionBegin;
410298fd71SBarry Smith   *xx = NULL;
428fe8eb27SJed Brown   if (!shell->left) {
438fe8eb27SJed Brown     *xx = x;
448fe8eb27SJed Brown   } else {
458fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
468fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
478fe8eb27SJed Brown     *xx  = shell->left_work;
488fe8eb27SJed Brown   }
498fe8eb27SJed Brown   PetscFunctionReturn(0);
508fe8eb27SJed Brown }
518fe8eb27SJed Brown 
520c0fd78eSBarry Smith /*
530c0fd78eSBarry Smith      xx = diag(right)*x
540c0fd78eSBarry Smith */
558fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
568fe8eb27SJed Brown {
578fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
588fe8eb27SJed Brown   PetscErrorCode ierr;
598fe8eb27SJed Brown 
608fe8eb27SJed Brown   PetscFunctionBegin;
610298fd71SBarry Smith   *xx = NULL;
628fe8eb27SJed Brown   if (!shell->right) {
638fe8eb27SJed Brown     *xx = x;
648fe8eb27SJed Brown   } else {
658fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
668fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
678fe8eb27SJed Brown     *xx  = shell->right_work;
688fe8eb27SJed Brown   }
698fe8eb27SJed Brown   PetscFunctionReturn(0);
708fe8eb27SJed Brown }
718fe8eb27SJed Brown 
720c0fd78eSBarry Smith /*
730c0fd78eSBarry Smith     x = diag(left)*x
740c0fd78eSBarry Smith */
758fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
768fe8eb27SJed Brown {
778fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
788fe8eb27SJed Brown   PetscErrorCode ierr;
798fe8eb27SJed Brown 
808fe8eb27SJed Brown   PetscFunctionBegin;
818fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
828fe8eb27SJed Brown   PetscFunctionReturn(0);
838fe8eb27SJed Brown }
848fe8eb27SJed Brown 
850c0fd78eSBarry Smith /*
860c0fd78eSBarry Smith     x = diag(right)*x
870c0fd78eSBarry Smith */
888fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
898fe8eb27SJed Brown {
908fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
918fe8eb27SJed Brown   PetscErrorCode ierr;
928fe8eb27SJed Brown 
938fe8eb27SJed Brown   PetscFunctionBegin;
948fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
958fe8eb27SJed Brown   PetscFunctionReturn(0);
968fe8eb27SJed Brown }
978fe8eb27SJed Brown 
980c0fd78eSBarry Smith /*
990c0fd78eSBarry Smith          Y = vscale*Y + diag(dshift)*X + vshift*X
1000c0fd78eSBarry Smith 
1010c0fd78eSBarry Smith          On input Y already contains A*x
1020c0fd78eSBarry Smith */
1038fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
1048fe8eb27SJed Brown {
1058fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1068fe8eb27SJed Brown   PetscErrorCode ierr;
1078fe8eb27SJed Brown 
1088fe8eb27SJed Brown   PetscFunctionBegin;
1098fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
1108fe8eb27SJed Brown     PetscInt          i,m;
1118fe8eb27SJed Brown     const PetscScalar *x,*d;
1128fe8eb27SJed Brown     PetscScalar       *y;
1138fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
1148fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1158fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1168fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1178fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
1188fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1198fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1208fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
1210c0fd78eSBarry Smith   } else {
122027c5fb5SJed Brown     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
1238fe8eb27SJed Brown   }
124d4c7de66SBarry Smith   if (shell->vshift != 0.0) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */
1258fe8eb27SJed Brown   PetscFunctionReturn(0);
1268fe8eb27SJed Brown }
127e51e0e81SBarry Smith 
1289d225801SJed Brown /*@
129a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
130b4fd4287SBarry Smith 
131c7fcc2eaSBarry Smith     Not Collective
132c7fcc2eaSBarry Smith 
133b4fd4287SBarry Smith     Input Parameter:
134b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
135b4fd4287SBarry Smith 
136b4fd4287SBarry Smith     Output Parameter:
137b4fd4287SBarry Smith .   ctx - the user provided context
138b4fd4287SBarry Smith 
13915091d37SBarry Smith     Level: advanced
14015091d37SBarry Smith 
14195452b02SPatrick Sanan    Fortran Notes:
14295452b02SPatrick Sanan     To use this from Fortran you must write a Fortran interface definition for this
143daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
144a62d957aSLois Curfman McInnes 
145a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context
146a62d957aSLois Curfman McInnes 
147ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
1480bc0a719SBarry Smith @*/
1498fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
150b4fd4287SBarry Smith {
151dfbe8321SBarry Smith   PetscErrorCode ierr;
152ace3abfcSBarry Smith   PetscBool      flg;
153273d9f13SBarry Smith 
1543a40ed3dSBarry Smith   PetscFunctionBegin;
1550700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1564482741eSBarry Smith   PetscValidPointer(ctx,2);
157251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
158940183f3SBarry Smith   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
159ce94432eSBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
1603a40ed3dSBarry Smith   PetscFunctionReturn(0);
161b4fd4287SBarry Smith }
162b4fd4287SBarry Smith 
163dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
164e51e0e81SBarry Smith {
165dfbe8321SBarry Smith   PetscErrorCode ierr;
166bf0cc555SLisandro Dalcin   Mat_Shell      *shell = (Mat_Shell*)mat->data;
167ed3cc1f0SBarry Smith 
1683a40ed3dSBarry Smith   PetscFunctionBegin;
16928f357e3SAlex Fikl   if (shell->ops->destroy) {
17028f357e3SAlex Fikl     ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr);
171bf0cc555SLisandro Dalcin   }
1720c0fd78eSBarry Smith   ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
1730c0fd78eSBarry Smith   ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
1740c0fd78eSBarry Smith   ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
1758fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
1768fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
1775edf6516SJed Brown   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
1785edf6516SJed Brown   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
1799f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
180bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
1813a40ed3dSBarry Smith   PetscFunctionReturn(0);
182e51e0e81SBarry Smith }
183e51e0e81SBarry Smith 
1847fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
1857fabda3fSAlex Fikl {
18628f357e3SAlex Fikl   Mat_Shell       *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
1877fabda3fSAlex Fikl   PetscErrorCode  ierr;
188cb8c8a77SBarry Smith   PetscBool       matflg;
1897fabda3fSAlex Fikl 
1907fabda3fSAlex Fikl   PetscFunctionBegin;
19128f357e3SAlex Fikl   ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr);
192cb8c8a77SBarry Smith   if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL");
1937fabda3fSAlex Fikl 
194cb8c8a77SBarry Smith   ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
195cb8c8a77SBarry Smith   ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
1967fabda3fSAlex Fikl 
197cb8c8a77SBarry Smith   if (shellA->ops->copy) {
19828f357e3SAlex Fikl     ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr);
199cb8c8a77SBarry Smith   }
2007fabda3fSAlex Fikl   shellB->vscale = shellA->vscale;
2017fabda3fSAlex Fikl   shellB->vshift = shellA->vshift;
2020c0fd78eSBarry Smith   if (shellA->dshift) {
2030c0fd78eSBarry Smith     if (!shellB->dshift) {
2040c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr);
2057fabda3fSAlex Fikl     }
2060c0fd78eSBarry Smith     ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr);
2077fabda3fSAlex Fikl   } else {
2089f137db4SBarry Smith     ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr);
2097fabda3fSAlex Fikl   }
2100c0fd78eSBarry Smith   if (shellA->left) {
2110c0fd78eSBarry Smith     if (!shellB->left) {
2120c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr);
2137fabda3fSAlex Fikl     }
2140c0fd78eSBarry Smith     ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr);
2157fabda3fSAlex Fikl   } else {
2169f137db4SBarry Smith     ierr = VecDestroy(&shellB->left);CHKERRQ(ierr);
2177fabda3fSAlex Fikl   }
2180c0fd78eSBarry Smith   if (shellA->right) {
2190c0fd78eSBarry Smith     if (!shellB->right) {
2200c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr);
2217fabda3fSAlex Fikl     }
2220c0fd78eSBarry Smith     ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr);
2237fabda3fSAlex Fikl   } else {
2249f137db4SBarry Smith     ierr = VecDestroy(&shellB->right);CHKERRQ(ierr);
2257fabda3fSAlex Fikl   }
22640e381c3SBarry Smith   ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr);
22740e381c3SBarry Smith   if (shellA->axpy) {
22840e381c3SBarry Smith     ierr                 = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr);
22940e381c3SBarry Smith     shellB->axpy        = shellA->axpy;
23040e381c3SBarry Smith     shellB->axpy_vscale = shellA->axpy_vscale;
23140e381c3SBarry Smith   }
2327fabda3fSAlex Fikl   PetscFunctionReturn(0);
2337fabda3fSAlex Fikl }
2347fabda3fSAlex Fikl 
235cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M)
236cb8c8a77SBarry Smith {
237cb8c8a77SBarry Smith   PetscErrorCode ierr;
238cb8c8a77SBarry Smith   void           *ctx;
239cb8c8a77SBarry Smith 
240cb8c8a77SBarry Smith   PetscFunctionBegin;
241cb8c8a77SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
242cb8c8a77SBarry Smith   ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr);
243cb8c8a77SBarry Smith   ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
244cb8c8a77SBarry Smith   PetscFunctionReturn(0);
245cb8c8a77SBarry Smith }
246cb8c8a77SBarry Smith 
247dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
248ef66eb69SBarry Smith {
249ef66eb69SBarry Smith   Mat_Shell        *shell = (Mat_Shell*)A->data;
250dfbe8321SBarry Smith   PetscErrorCode   ierr;
25125578ef6SJed Brown   Vec              xx;
252e3f487b0SBarry Smith   PetscObjectState instate,outstate;
253ef66eb69SBarry Smith 
254ef66eb69SBarry Smith   PetscFunctionBegin;
2558fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
256e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
257976e8c5aSLisandro Dalcin   if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL");
25828f357e3SAlex Fikl   ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr);
259e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
260e3f487b0SBarry Smith   if (instate == outstate) {
261e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
262e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
263e3f487b0SBarry Smith   }
2648fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2658fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
2669f137db4SBarry Smith 
2679f137db4SBarry Smith   if (shell->axpy) {
2689f137db4SBarry Smith     if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);}
2699f137db4SBarry Smith     ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr);
2709f137db4SBarry Smith     ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
2719f137db4SBarry Smith   }
272ef66eb69SBarry Smith   PetscFunctionReturn(0);
273ef66eb69SBarry Smith }
274ef66eb69SBarry Smith 
2755edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
2765edf6516SJed Brown {
2775edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
2785edf6516SJed Brown   PetscErrorCode ierr;
2795edf6516SJed Brown 
2805edf6516SJed Brown   PetscFunctionBegin;
2815edf6516SJed Brown   if (y == z) {
2825edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
2835edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
284b8430d49SBarry Smith     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
2855edf6516SJed Brown   } else {
2865edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
2875edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
2885edf6516SJed Brown   }
2895edf6516SJed Brown   PetscFunctionReturn(0);
2905edf6516SJed Brown }
2915edf6516SJed Brown 
29218be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
29318be62a5SSatish Balay {
29418be62a5SSatish Balay   Mat_Shell        *shell = (Mat_Shell*)A->data;
29518be62a5SSatish Balay   PetscErrorCode   ierr;
29625578ef6SJed Brown   Vec              xx;
297e3f487b0SBarry Smith   PetscObjectState instate,outstate;
29818be62a5SSatish Balay 
29918be62a5SSatish Balay   PetscFunctionBegin;
3008fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
301e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
3020c0fd78eSBarry Smith   if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL");
30328f357e3SAlex Fikl   ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr);
304e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
305e3f487b0SBarry Smith   if (instate == outstate) {
306e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
307e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
308e3f487b0SBarry Smith   }
3098fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
3108fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
31118be62a5SSatish Balay   PetscFunctionReturn(0);
31218be62a5SSatish Balay }
31318be62a5SSatish Balay 
3145edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3155edf6516SJed Brown {
3165edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3175edf6516SJed Brown   PetscErrorCode ierr;
3185edf6516SJed Brown 
3195edf6516SJed Brown   PetscFunctionBegin;
3205edf6516SJed Brown   if (y == z) {
3215edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
3225edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
323*b253ae0bS“Marek     ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr);
3245edf6516SJed Brown   } else {
3255edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
3265edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3275edf6516SJed Brown   }
3285edf6516SJed Brown   PetscFunctionReturn(0);
3295edf6516SJed Brown }
3305edf6516SJed Brown 
3310c0fd78eSBarry Smith /*
3320c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
3330c0fd78eSBarry Smith */
33418be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
33518be62a5SSatish Balay {
33618be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
33718be62a5SSatish Balay   PetscErrorCode ierr;
33818be62a5SSatish Balay 
33918be62a5SSatish Balay   PetscFunctionBegin;
3400c0fd78eSBarry Smith   if (shell->ops->getdiagonal) {
34128f357e3SAlex Fikl     ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr);
342305211d5SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Must provide shell matrix with routine to return diagonal using\nMatShellSetOperation(S,MATOP_GET_DIAGONAL,...)");
34318be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
3448fe8eb27SJed Brown   if (shell->dshift) {
3450c0fd78eSBarry Smith     ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr);
34618be62a5SSatish Balay   }
3470c0fd78eSBarry Smith   ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
3488fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
3498fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
35081c02519SBarry Smith   if (shell->axpy) {
35181c02519SBarry Smith     if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);}
35281c02519SBarry Smith     ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr);
35381c02519SBarry Smith     ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
35481c02519SBarry Smith   }
35518be62a5SSatish Balay   PetscFunctionReturn(0);
35618be62a5SSatish Balay }
35718be62a5SSatish Balay 
358f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
359ef66eb69SBarry Smith {
360ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3618fe8eb27SJed Brown   PetscErrorCode ierr;
362b24ad042SBarry Smith 
363ef66eb69SBarry Smith   PetscFunctionBegin;
3640c0fd78eSBarry Smith   if (shell->left || shell->right) {
3658fe8eb27SJed Brown     if (!shell->dshift) {
3660c0fd78eSBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
3670c0fd78eSBarry Smith       ierr = VecSet(shell->dshift,a);CHKERRQ(ierr);
3680c0fd78eSBarry Smith     } else {
3690c0fd78eSBarry Smith       if (shell->left)  {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
3700c0fd78eSBarry Smith       if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
3719f137db4SBarry Smith       ierr = VecShift(shell->dshift,a);CHKERRQ(ierr);
3720c0fd78eSBarry Smith     }
3738fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
3748fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
3758fe8eb27SJed Brown   } else shell->vshift += a;
376ef66eb69SBarry Smith   PetscFunctionReturn(0);
377ef66eb69SBarry Smith }
378ef66eb69SBarry Smith 
379*b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s)
3806464bf51SAlex Fikl {
3816464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
3826464bf51SAlex Fikl   PetscErrorCode ierr;
3836464bf51SAlex Fikl 
3846464bf51SAlex Fikl   PetscFunctionBegin;
3850c0fd78eSBarry Smith   if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);}
3860c0fd78eSBarry Smith   if (shell->left || shell->right) {
38792fabd1fSBarry Smith     if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);}
3889f137db4SBarry Smith     if (shell->left && shell->right)  {
3899f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
3909f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
3919f137db4SBarry Smith     } else if (shell->left) {
3929f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
3939f137db4SBarry Smith     } else {
3949f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr);
3959f137db4SBarry Smith     }
3969f137db4SBarry Smith     if (!shell->dshift) {
3979f137db4SBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
3989f137db4SBarry Smith       ierr = VecCopy(shell->dshift,shell->right_work);CHKERRQ(ierr);
3999f137db4SBarry Smith     } else {
400*b253ae0bS“Marek       ierr = VecAXPY(shell->dshift,s,shell->right_work);CHKERRQ(ierr);
4019f137db4SBarry Smith     }
4020c0fd78eSBarry Smith   } else {
403*b253ae0bS“Marek     ierr = VecAXPY(shell->dshift,s,D);CHKERRQ(ierr);
404*b253ae0bS“Marek   }
405*b253ae0bS“Marek   PetscFunctionReturn(0);
406*b253ae0bS“Marek }
407*b253ae0bS“Marek 
408*b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
409*b253ae0bS“Marek {
410*b253ae0bS“Marek   Vec            d;
411*b253ae0bS“Marek   PetscErrorCode ierr;
412*b253ae0bS“Marek 
413*b253ae0bS“Marek   PetscFunctionBegin;
414*b253ae0bS“Marek   if (ins == INSERT_VALUES) {
415*b253ae0bS“Marek     if (!A->ops->getdiagonal) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first");
416*b253ae0bS“Marek     ierr = VecDuplicate(D,&d);CHKERRQ(ierr);
417*b253ae0bS“Marek     ierr = MatGetDiagonal(A,d);CHKERRQ(ierr);
418*b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,d,-1.);CHKERRQ(ierr);
419*b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr);
420*b253ae0bS“Marek     ierr = VecDestroy(&d);CHKERRQ(ierr);
421*b253ae0bS“Marek   } else {
422*b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr);
4236464bf51SAlex Fikl   }
4246464bf51SAlex Fikl   PetscFunctionReturn(0);
4256464bf51SAlex Fikl }
4266464bf51SAlex Fikl 
427f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
428ef66eb69SBarry Smith {
429ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4308fe8eb27SJed Brown   PetscErrorCode ierr;
431b24ad042SBarry Smith 
432ef66eb69SBarry Smith   PetscFunctionBegin;
433f4df32b1SMatthew Knepley   shell->vscale *= a;
4340c0fd78eSBarry Smith   shell->vshift *= a;
4352205254eSKarl Rupp   if (shell->dshift) {
4362205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
4370c0fd78eSBarry Smith   }
43881c02519SBarry Smith   shell->axpy_vscale *= a;
4398fe8eb27SJed Brown   PetscFunctionReturn(0);
44018be62a5SSatish Balay }
4418fe8eb27SJed Brown 
4428fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
4438fe8eb27SJed Brown {
4448fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4458fe8eb27SJed Brown   PetscErrorCode ierr;
4468fe8eb27SJed Brown 
4478fe8eb27SJed Brown   PetscFunctionBegin;
44881c02519SBarry Smith   if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation");
4498fe8eb27SJed Brown   if (left) {
4500c0fd78eSBarry Smith     if (!shell->left) {
4510c0fd78eSBarry Smith       ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr);
4528fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
4530c0fd78eSBarry Smith     } else {
4540c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
45518be62a5SSatish Balay     }
456ef66eb69SBarry Smith   }
4578fe8eb27SJed Brown   if (right) {
4580c0fd78eSBarry Smith     if (!shell->right) {
4590c0fd78eSBarry Smith       ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr);
4608fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
4610c0fd78eSBarry Smith     } else {
4620c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
4638fe8eb27SJed Brown     }
4648fe8eb27SJed Brown   }
465ef66eb69SBarry Smith   PetscFunctionReturn(0);
466ef66eb69SBarry Smith }
467ef66eb69SBarry Smith 
468dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
469ef66eb69SBarry Smith {
470ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4710c0fd78eSBarry Smith   PetscErrorCode ierr;
472ef66eb69SBarry Smith 
473ef66eb69SBarry Smith   PetscFunctionBegin;
4748fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
475ef66eb69SBarry Smith     shell->vshift = 0.0;
476ef66eb69SBarry Smith     shell->vscale = 1.0;
4770c0fd78eSBarry Smith     ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
4780c0fd78eSBarry Smith     ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
4790c0fd78eSBarry Smith     ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
48081c02519SBarry Smith     ierr = MatDestroy(&shell->axpy);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,
52400849d43SBarry Smith                                        0,
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 
71595452b02SPatrick Sanan    Fortran Notes:
71695452b02SPatrick Sanan     To use this from Fortran with a ctx you must write an interface definition for this
717daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
718daf670e6SBarry Smith     in as the ctx argument.
719f38a8d56SBarry Smith 
720f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
721f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
722645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
723645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
724645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
725645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
726645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
727645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
728645985a0SLois Curfman McInnes    For example,
729f39d1f56SLois Curfman McInnes 
730f39d1f56SLois Curfman McInnes $
731f39d1f56SLois Curfman McInnes $     Vec x, y
7327b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
733645985a0SLois Curfman McInnes $     Mat A
734f39d1f56SLois Curfman McInnes $
735c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
736c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
737f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
738c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
739c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
740c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
741645985a0SLois Curfman McInnes $     MatMult(A,x,y);
742645985a0SLois Curfman McInnes $     MatDestroy(A);
743f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
744645985a0SLois Curfman McInnes $
745e51e0e81SBarry Smith 
7469b53d762SBarry Smith 
7479b53d762SBarry Smith    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), and MatScale() internally so these
7489b53d762SBarry Smith    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.
7499b53d762SBarry Smith 
7509b53d762SBarry Smith 
7510c0fd78eSBarry Smith     For rectangular matrices do all the scalings and shifts make sense?
7520c0fd78eSBarry Smith 
75395452b02SPatrick Sanan     Developers Notes:
75495452b02SPatrick Sanan     Regarding shifting and scaling. The general form is
7550c0fd78eSBarry Smith 
7560c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
7570c0fd78eSBarry Smith 
7580c0fd78eSBarry Smith       The order you apply the operations is important. For example if you have a dshift then
7590c0fd78eSBarry Smith       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
7600c0fd78eSBarry Smith       you get s*vscale*A + diag(shift)
7610c0fd78eSBarry Smith 
7620c0fd78eSBarry Smith           A is the user provided function.
7630c0fd78eSBarry Smith 
7640b627109SLois Curfman McInnes .keywords: matrix, shell, create
7650b627109SLois Curfman McInnes 
7660c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
767e51e0e81SBarry Smith @*/
7687087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
769e51e0e81SBarry Smith {
770dfbe8321SBarry Smith   PetscErrorCode ierr;
771ed3cc1f0SBarry Smith 
7723a40ed3dSBarry Smith   PetscFunctionBegin;
773f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
774f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
775273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
776273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
7770fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
778273d9f13SBarry Smith   PetscFunctionReturn(0);
779c7fcc2eaSBarry Smith }
780c7fcc2eaSBarry Smith 
781c6866cfdSSatish Balay /*@
782273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
783c7fcc2eaSBarry Smith 
7843f9fe445SBarry Smith    Logically Collective on Mat
785c7fcc2eaSBarry Smith 
786273d9f13SBarry Smith     Input Parameters:
787273d9f13SBarry Smith +   mat - the shell matrix
788273d9f13SBarry Smith -   ctx - the context
789273d9f13SBarry Smith 
790273d9f13SBarry Smith    Level: advanced
791273d9f13SBarry Smith 
79295452b02SPatrick Sanan    Fortran Notes:
79395452b02SPatrick Sanan     To use this from Fortran you must write a Fortran interface definition for this
794daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
795273d9f13SBarry Smith 
796273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
7970bc0a719SBarry Smith @*/
7987087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
799273d9f13SBarry Smith {
800273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
801dfbe8321SBarry Smith   PetscErrorCode ierr;
802ace3abfcSBarry Smith   PetscBool      flg;
803273d9f13SBarry Smith 
804273d9f13SBarry Smith   PetscFunctionBegin;
8050700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
806251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
807273d9f13SBarry Smith   if (flg) {
808273d9f13SBarry Smith     shell->ctx = ctx;
809ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
8103a40ed3dSBarry Smith   PetscFunctionReturn(0);
811e51e0e81SBarry Smith }
812e51e0e81SBarry Smith 
8130c0fd78eSBarry Smith /*@
8140c0fd78eSBarry Smith     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
8150c0fd78eSBarry Smith           after MatCreateShell()
8160c0fd78eSBarry Smith 
8170c0fd78eSBarry Smith    Logically Collective on Mat
8180c0fd78eSBarry Smith 
8190c0fd78eSBarry Smith     Input Parameter:
8200c0fd78eSBarry Smith .   mat - the shell matrix
8210c0fd78eSBarry Smith 
8220c0fd78eSBarry Smith   Level: advanced
8230c0fd78eSBarry Smith 
8240c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
8250c0fd78eSBarry Smith @*/
8260c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A)
8270c0fd78eSBarry Smith {
8280c0fd78eSBarry Smith   PetscErrorCode ierr;
829976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
8300c0fd78eSBarry Smith   PetscBool      flg;
8310c0fd78eSBarry Smith 
8320c0fd78eSBarry Smith   PetscFunctionBegin;
8330c0fd78eSBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
8340c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr);
8350c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
836976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)A->data;
8370c0fd78eSBarry Smith   shell->managescalingshifts = PETSC_FALSE;
838976e8c5aSLisandro Dalcin   A->ops->diagonalset   = NULL;
839976e8c5aSLisandro Dalcin   A->ops->diagonalscale = NULL;
840976e8c5aSLisandro Dalcin   A->ops->scale         = NULL;
841976e8c5aSLisandro Dalcin   A->ops->shift         = NULL;
842976e8c5aSLisandro Dalcin   A->ops->axpy          = NULL;
8430c0fd78eSBarry Smith   PetscFunctionReturn(0);
8440c0fd78eSBarry Smith }
8450c0fd78eSBarry Smith 
846c16cb8f2SBarry Smith /*@C
847f3b1f45cSBarry Smith     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
848f3b1f45cSBarry Smith 
849f3b1f45cSBarry Smith    Logically Collective on Mat
850f3b1f45cSBarry Smith 
851f3b1f45cSBarry Smith     Input Parameters:
852f3b1f45cSBarry Smith +   mat - the shell matrix
853f3b1f45cSBarry Smith .   f - the function
854f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
855f3b1f45cSBarry Smith -   ctx - an optional context for the function
856f3b1f45cSBarry Smith 
857f3b1f45cSBarry Smith    Output Parameter:
858f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
859f3b1f45cSBarry Smith 
860f3b1f45cSBarry Smith    Options Database:
861f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
862f3b1f45cSBarry Smith 
863f3b1f45cSBarry Smith    Level: advanced
864f3b1f45cSBarry Smith 
86595452b02SPatrick Sanan    Fortran Notes:
86695452b02SPatrick Sanan     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) {
896fc7aafd1SBarry 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));CHKERRQ(ierr);
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) {
902fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
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 
93195452b02SPatrick Sanan    Fortran Notes:
93295452b02SPatrick Sanan     Not supported from Fortran
933f3b1f45cSBarry Smith 
934f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
935f3b1f45cSBarry Smith @*/
936f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
937f3b1f45cSBarry Smith {
938f3b1f45cSBarry Smith   PetscErrorCode ierr;
939f3b1f45cSBarry Smith   Vec            x,y,z;
940f3b1f45cSBarry Smith   PetscInt       m,n,M,N;
941f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
942f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
94374e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
944f3b1f45cSBarry Smith 
945f3b1f45cSBarry Smith   PetscFunctionBegin;
946f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
947f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
948f3b1f45cSBarry Smith   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
949f3b1f45cSBarry Smith   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
950f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
951f3b1f45cSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
952f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
953f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
954f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
955f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
956f3b1f45cSBarry Smith   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
957f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr);
958f3b1f45cSBarry Smith 
959f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
960f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
961f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
962f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
963f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
964f3b1f45cSBarry Smith     flag = PETSC_FALSE;
965f3b1f45cSBarry Smith     if (v) {
966fc7aafd1SBarry 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));CHKERRQ(ierr);
967f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
968f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
969f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
970f3b1f45cSBarry Smith     }
971f3b1f45cSBarry Smith   } else if (v) {
972fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
973f3b1f45cSBarry Smith   }
974f3b1f45cSBarry Smith   if (flg) *flg = flag;
975f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
976f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
977f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
978f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
979f3b1f45cSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
980f3b1f45cSBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
981f3b1f45cSBarry Smith   ierr = VecDestroy(&z);CHKERRQ(ierr);
982f3b1f45cSBarry Smith   PetscFunctionReturn(0);
983f3b1f45cSBarry Smith }
984f3b1f45cSBarry Smith 
985f3b1f45cSBarry Smith /*@C
9860c0fd78eSBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
987e51e0e81SBarry Smith 
9883f9fe445SBarry Smith    Logically Collective on Mat
989fee21e36SBarry Smith 
990c7fcc2eaSBarry Smith     Input Parameters:
991c7fcc2eaSBarry Smith +   mat - the shell matrix
992c7fcc2eaSBarry Smith .   op - the name of the operation
993c7fcc2eaSBarry Smith -   f - the function that provides the operation.
994c7fcc2eaSBarry Smith 
99515091d37SBarry Smith    Level: advanced
99615091d37SBarry Smith 
997fae171e0SBarry Smith     Usage:
998e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
999f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
1000c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
10010b627109SLois Curfman McInnes 
1002a62d957aSLois Curfman McInnes     Notes:
1003e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
10041c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
1005a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
10061c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
1007a62d957aSLois Curfman McInnes 
1008a4d85fd7SBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
1009deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
1010deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
1011deebb3c3SLois Curfman McInnes     routines, e.g.,
1012a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1013a62d957aSLois Curfman McInnes 
10144aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
10154aa34b0aSBarry Smith     nonzero on failure.
10164aa34b0aSBarry Smith 
1017a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
1018a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
1019a62d957aSLois Curfman McInnes     set by MatCreateShell().
1020a62d957aSLois Curfman McInnes 
102195452b02SPatrick Sanan     Fortran Notes:
102295452b02SPatrick Sanan     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 {
1033ace3abfcSBarry Smith   PetscBool      flg;
1034976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
1035976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
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");
1041976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)mat->data;
1042976e8c5aSLisandro Dalcin 
10435edf6516SJed Brown   switch (op) {
10445edf6516SJed Brown   case MATOP_DESTROY:
104528f357e3SAlex Fikl     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
10465edf6516SJed Brown     break;
1047976e8c5aSLisandro Dalcin   case MATOP_VIEW:
1048976e8c5aSLisandro Dalcin     if (!mat->ops->viewnative) {
1049976e8c5aSLisandro Dalcin       mat->ops->viewnative = mat->ops->view;
1050976e8c5aSLisandro Dalcin     }
1051976e8c5aSLisandro Dalcin     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
1052976e8c5aSLisandro Dalcin     break;
1053976e8c5aSLisandro Dalcin   case MATOP_COPY:
1054976e8c5aSLisandro Dalcin     shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f;
1055976e8c5aSLisandro Dalcin     break;
10566464bf51SAlex Fikl   case MATOP_DIAGONAL_SET:
10570c0fd78eSBarry Smith   case MATOP_DIAGONAL_SCALE:
10580c0fd78eSBarry Smith   case MATOP_SHIFT:
10590c0fd78eSBarry Smith   case MATOP_SCALE:
10609f137db4SBarry Smith   case MATOP_AXPY:
10610c0fd78eSBarry Smith     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
10620c0fd78eSBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
10636464bf51SAlex Fikl     break;
10640c0fd78eSBarry Smith   case MATOP_GET_DIAGONAL:
1065976e8c5aSLisandro Dalcin     if (shell->managescalingshifts) {
1066976e8c5aSLisandro Dalcin       shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1067976e8c5aSLisandro Dalcin       mat->ops->getdiagonal   = MatGetDiagonal_Shell;
1068976e8c5aSLisandro Dalcin     } else {
1069976e8c5aSLisandro Dalcin       shell->ops->getdiagonal = NULL;
1070976e8c5aSLisandro Dalcin       mat->ops->getdiagonal   = (PetscErrorCode (*)(Mat,Vec))f;
107140a5a6c4SBarry Smith     }
10725edf6516SJed Brown     break;
10735edf6516SJed Brown   case MATOP_MULT:
10749f137db4SBarry Smith     if (shell->managescalingshifts) {
10759f137db4SBarry Smith       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10769f137db4SBarry Smith       mat->ops->mult   = MatMult_Shell;
1077976e8c5aSLisandro Dalcin     } else {
1078976e8c5aSLisandro Dalcin       shell->ops->mult = NULL;
1079976e8c5aSLisandro Dalcin       mat->ops->mult   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1080976e8c5aSLisandro Dalcin     }
10815edf6516SJed Brown     break;
10825edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
10839f137db4SBarry Smith     if (shell->managescalingshifts) {
10849f137db4SBarry Smith       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10855259c5a4SBarry Smith       mat->ops->multtranspose   = MatMultTranspose_Shell;
1086976e8c5aSLisandro Dalcin     } else {
1087976e8c5aSLisandro Dalcin       shell->ops->multtranspose = NULL;
1088976e8c5aSLisandro Dalcin       mat->ops->multtranspose   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1089976e8c5aSLisandro Dalcin     }
10905edf6516SJed Brown     break;
10915edf6516SJed Brown   default:
10925edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
10935ab264f3SAlp Dener     break;
1094a62d957aSLois Curfman McInnes   }
10953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1096e51e0e81SBarry Smith }
1097f0479e8cSBarry Smith 
1098d4bb536fSBarry Smith /*@C
1099d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
1100d4bb536fSBarry Smith 
1101c7fcc2eaSBarry Smith     Not Collective
1102c7fcc2eaSBarry Smith 
1103d4bb536fSBarry Smith     Input Parameters:
1104c7fcc2eaSBarry Smith +   mat - the shell matrix
1105c7fcc2eaSBarry Smith -   op - the name of the operation
1106d4bb536fSBarry Smith 
1107d4bb536fSBarry Smith     Output Parameter:
1108d4bb536fSBarry Smith .   f - the function that provides the operation.
1109d4bb536fSBarry Smith 
111015091d37SBarry Smith     Level: advanced
111115091d37SBarry Smith 
1112d4bb536fSBarry Smith     Notes:
1113e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
1114d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
1115d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
1116d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
1117d4bb536fSBarry Smith 
1118d4bb536fSBarry Smith     All user-provided functions have the same calling
1119d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
1120d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
1121d4bb536fSBarry Smith     routines, e.g.,
1122d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1123d4bb536fSBarry Smith 
1124d4bb536fSBarry Smith     Within each user-defined routine, the user should call
1125d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
1126d4bb536fSBarry Smith     set by MatCreateShell().
1127d4bb536fSBarry Smith 
1128d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
1129d4bb536fSBarry Smith 
1130ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1131d4bb536fSBarry Smith @*/
11327087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
1133d4bb536fSBarry Smith {
1134ace3abfcSBarry Smith   PetscBool      flg;
1135976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
1136976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
1137273d9f13SBarry Smith 
11383a40ed3dSBarry Smith   PetscFunctionBegin;
11390700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1140976e8c5aSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1141976e8c5aSLisandro Dalcin   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1142976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)mat->data;
1143976e8c5aSLisandro Dalcin 
114428f357e3SAlex Fikl   switch (op) {
114528f357e3SAlex Fikl   case MATOP_DESTROY:
114628f357e3SAlex Fikl     *f = (void (*)(void))shell->ops->destroy;
114728f357e3SAlex Fikl     break;
114828f357e3SAlex Fikl   case MATOP_VIEW:
114928f357e3SAlex Fikl     *f = (void (*)(void))mat->ops->view;
115028f357e3SAlex Fikl     break;
1151976e8c5aSLisandro Dalcin   case MATOP_COPY:
1152976e8c5aSLisandro Dalcin     *f = (void (*)(void))shell->ops->copy;
1153976e8c5aSLisandro Dalcin     break;
1154976e8c5aSLisandro Dalcin   case MATOP_DIAGONAL_SET:
1155976e8c5aSLisandro Dalcin   case MATOP_DIAGONAL_SCALE:
1156976e8c5aSLisandro Dalcin   case MATOP_SHIFT:
1157976e8c5aSLisandro Dalcin   case MATOP_SCALE:
1158976e8c5aSLisandro Dalcin   case MATOP_AXPY:
1159976e8c5aSLisandro Dalcin     *f = (((void (**)(void))mat->ops)[op]);
1160976e8c5aSLisandro Dalcin     break;
1161976e8c5aSLisandro Dalcin   case MATOP_GET_DIAGONAL:
1162976e8c5aSLisandro Dalcin     if (shell->ops->getdiagonal)
1163976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->getdiagonal;
1164976e8c5aSLisandro Dalcin     else
1165976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1166976e8c5aSLisandro Dalcin     break;
1167976e8c5aSLisandro Dalcin   case MATOP_MULT:
1168976e8c5aSLisandro Dalcin     if (shell->ops->mult)
1169976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->mult;
1170976e8c5aSLisandro Dalcin     else
1171976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1172976e8c5aSLisandro Dalcin     break;
1173976e8c5aSLisandro Dalcin   case MATOP_MULT_TRANSPOSE:
1174976e8c5aSLisandro Dalcin     if (shell->ops->multtranspose)
1175976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->multtranspose;
1176976e8c5aSLisandro Dalcin     else
1177976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1178976e8c5aSLisandro Dalcin     break;
117928f357e3SAlex Fikl   default:
1180c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
1181d4bb536fSBarry Smith   }
11823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1183d4bb536fSBarry Smith }
1184