xref: /petsc/src/mat/impls/shell/shell.c (revision ad2f5c3f9718753398e340b8a0bf2bce5a4f6f72)
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);
323dac989eaS“Marek     ierr = VecAXPY(z,1.0,shell->left_add_work);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 
379b253ae0bS“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 {
400b253ae0bS“Marek       ierr = VecAXPY(shell->dshift,s,shell->right_work);CHKERRQ(ierr);
4019f137db4SBarry Smith     }
4020c0fd78eSBarry Smith   } else {
403b253ae0bS“Marek     ierr = VecAXPY(shell->dshift,s,D);CHKERRQ(ierr);
404b253ae0bS“Marek   }
405b253ae0bS“Marek   PetscFunctionReturn(0);
406b253ae0bS“Marek }
407b253ae0bS“Marek 
408b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
409b253ae0bS“Marek {
410b253ae0bS“Marek   Vec            d;
411b253ae0bS“Marek   PetscErrorCode ierr;
412b253ae0bS“Marek 
413b253ae0bS“Marek   PetscFunctionBegin;
414b253ae0bS“Marek   if (ins == INSERT_VALUES) {
415b253ae0bS“Marek     if (!A->ops->getdiagonal) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first");
416b253ae0bS“Marek     ierr = VecDuplicate(D,&d);CHKERRQ(ierr);
417b253ae0bS“Marek     ierr = MatGetDiagonal(A,d);CHKERRQ(ierr);
418b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,d,-1.);CHKERRQ(ierr);
419b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr);
420b253ae0bS“Marek     ierr = VecDestroy(&d);CHKERRQ(ierr);
421b253ae0bS“Marek   } else {
422b253ae0bS“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 
4853b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
4863b49f96aSBarry Smith {
4873b49f96aSBarry Smith   PetscFunctionBegin;
4883b49f96aSBarry Smith   *missing = PETSC_FALSE;
4893b49f96aSBarry Smith   PetscFunctionReturn(0);
4903b49f96aSBarry Smith }
4913b49f96aSBarry Smith 
4929f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
4939f137db4SBarry Smith {
4949f137db4SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4959f137db4SBarry Smith   PetscErrorCode ierr;
4969f137db4SBarry Smith 
4979f137db4SBarry Smith   PetscFunctionBegin;
4989f137db4SBarry Smith   ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr);
4999f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
5009f137db4SBarry Smith   shell->axpy        = X;
5019f137db4SBarry Smith   shell->axpy_vscale = a;
5029f137db4SBarry Smith   PetscFunctionReturn(0);
5039f137db4SBarry Smith }
5049f137db4SBarry Smith 
50509dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
50620563c6bSBarry Smith                                        0,
50720563c6bSBarry Smith                                        0,
5089f137db4SBarry Smith                                        0,
5090c0fd78eSBarry Smith                                 /* 4*/ MatMultAdd_Shell,
5109f137db4SBarry Smith                                        0,
5110c0fd78eSBarry Smith                                        MatMultTransposeAdd_Shell,
512b951964fSBarry Smith                                        0,
513b951964fSBarry Smith                                        0,
514b951964fSBarry Smith                                        0,
51597304618SKris Buschelman                                 /*10*/ 0,
516b951964fSBarry Smith                                        0,
517b951964fSBarry Smith                                        0,
518b951964fSBarry Smith                                        0,
519b951964fSBarry Smith                                        0,
52097304618SKris Buschelman                                 /*15*/ 0,
521b951964fSBarry Smith                                        0,
52200849d43SBarry Smith                                        0,
5238fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
524b951964fSBarry Smith                                        0,
52597304618SKris Buschelman                                 /*20*/ 0,
526ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
527b951964fSBarry Smith                                        0,
528b951964fSBarry Smith                                        0,
529d519adbfSMatthew Knepley                                 /*24*/ 0,
530b951964fSBarry Smith                                        0,
531b951964fSBarry Smith                                        0,
532b951964fSBarry Smith                                        0,
533b951964fSBarry Smith                                        0,
534d519adbfSMatthew Knepley                                 /*29*/ 0,
535b951964fSBarry Smith                                        0,
536273d9f13SBarry Smith                                        0,
537b951964fSBarry Smith                                        0,
538b951964fSBarry Smith                                        0,
539cb8c8a77SBarry Smith                                 /*34*/ MatDuplicate_Shell,
540b951964fSBarry Smith                                        0,
541b951964fSBarry Smith                                        0,
54209dc0095SBarry Smith                                        0,
54309dc0095SBarry Smith                                        0,
5449f137db4SBarry Smith                                 /*39*/ MatAXPY_Shell,
54509dc0095SBarry Smith                                        0,
54609dc0095SBarry Smith                                        0,
54709dc0095SBarry Smith                                        0,
548cb8c8a77SBarry Smith                                        MatCopy_Shell,
549d519adbfSMatthew Knepley                                 /*44*/ 0,
550ef66eb69SBarry Smith                                        MatScale_Shell,
551ef66eb69SBarry Smith                                        MatShift_Shell,
5526464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
55309dc0095SBarry Smith                                        0,
554f73d5cc4SBarry Smith                                 /*49*/ 0,
55509dc0095SBarry Smith                                        0,
55609dc0095SBarry Smith                                        0,
55709dc0095SBarry Smith                                        0,
55809dc0095SBarry Smith                                        0,
559d519adbfSMatthew Knepley                                 /*54*/ 0,
56009dc0095SBarry Smith                                        0,
56109dc0095SBarry Smith                                        0,
56209dc0095SBarry Smith                                        0,
56309dc0095SBarry Smith                                        0,
564d519adbfSMatthew Knepley                                 /*59*/ 0,
565b9b97703SBarry Smith                                        MatDestroy_Shell,
56609dc0095SBarry Smith                                        0,
567357abbc8SBarry Smith                                        0,
568273d9f13SBarry Smith                                        0,
569d519adbfSMatthew Knepley                                 /*64*/ 0,
570273d9f13SBarry Smith                                        0,
571273d9f13SBarry Smith                                        0,
572273d9f13SBarry Smith                                        0,
573273d9f13SBarry Smith                                        0,
574d519adbfSMatthew Knepley                                 /*69*/ 0,
575273d9f13SBarry Smith                                        0,
576c87e5d42SMatthew Knepley                                        MatConvert_Shell,
577273d9f13SBarry Smith                                        0,
57897304618SKris Buschelman                                        0,
579d519adbfSMatthew Knepley                                 /*74*/ 0,
58097304618SKris Buschelman                                        0,
58197304618SKris Buschelman                                        0,
58297304618SKris Buschelman                                        0,
58397304618SKris Buschelman                                        0,
584d519adbfSMatthew Knepley                                 /*79*/ 0,
58597304618SKris Buschelman                                        0,
58697304618SKris Buschelman                                        0,
58797304618SKris Buschelman                                        0,
588865e5f61SKris Buschelman                                        0,
589d519adbfSMatthew Knepley                                 /*84*/ 0,
590865e5f61SKris Buschelman                                        0,
591865e5f61SKris Buschelman                                        0,
592865e5f61SKris Buschelman                                        0,
593865e5f61SKris Buschelman                                        0,
594d519adbfSMatthew Knepley                                 /*89*/ 0,
595865e5f61SKris Buschelman                                        0,
596865e5f61SKris Buschelman                                        0,
597865e5f61SKris Buschelman                                        0,
598865e5f61SKris Buschelman                                        0,
599d519adbfSMatthew Knepley                                 /*94*/ 0,
600865e5f61SKris Buschelman                                        0,
601865e5f61SKris Buschelman                                        0,
6023964eb88SJed Brown                                        0,
6033964eb88SJed Brown                                        0,
6043964eb88SJed Brown                                 /*99*/ 0,
6053964eb88SJed Brown                                        0,
6063964eb88SJed Brown                                        0,
6073964eb88SJed Brown                                        0,
6083964eb88SJed Brown                                        0,
6093964eb88SJed Brown                                /*104*/ 0,
6103964eb88SJed Brown                                        0,
6113964eb88SJed Brown                                        0,
6123964eb88SJed Brown                                        0,
6133964eb88SJed Brown                                        0,
6143964eb88SJed Brown                                /*109*/ 0,
6153964eb88SJed Brown                                        0,
6163964eb88SJed Brown                                        0,
6173964eb88SJed Brown                                        0,
6183b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
6193964eb88SJed Brown                                /*114*/ 0,
6203964eb88SJed Brown                                        0,
6213964eb88SJed Brown                                        0,
6223964eb88SJed Brown                                        0,
6233964eb88SJed Brown                                        0,
6243964eb88SJed Brown                                /*119*/ 0,
6253964eb88SJed Brown                                        0,
6263964eb88SJed Brown                                        0,
6273964eb88SJed Brown                                        0,
6283964eb88SJed Brown                                        0,
6293964eb88SJed Brown                                /*124*/ 0,
6303964eb88SJed Brown                                        0,
6313964eb88SJed Brown                                        0,
6323964eb88SJed Brown                                        0,
6333964eb88SJed Brown                                        0,
6343964eb88SJed Brown                                /*129*/ 0,
6353964eb88SJed Brown                                        0,
6363964eb88SJed Brown                                        0,
6373964eb88SJed Brown                                        0,
6383964eb88SJed Brown                                        0,
6393964eb88SJed Brown                                /*134*/ 0,
6403964eb88SJed Brown                                        0,
6413964eb88SJed Brown                                        0,
6423964eb88SJed Brown                                        0,
6433964eb88SJed Brown                                        0,
6443964eb88SJed Brown                                /*139*/ 0,
645f9426fe0SMark Adams                                        0,
6463964eb88SJed Brown                                        0
6473964eb88SJed Brown };
648273d9f13SBarry Smith 
6490bad9183SKris Buschelman /*MC
650fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
6510bad9183SKris Buschelman 
6520bad9183SKris Buschelman   Level: advanced
6530bad9183SKris Buschelman 
6540c0fd78eSBarry Smith .seealso: MatCreateShell()
6550bad9183SKris Buschelman M*/
6560bad9183SKris Buschelman 
6578cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
658273d9f13SBarry Smith {
659273d9f13SBarry Smith   Mat_Shell      *b;
660dfbe8321SBarry Smith   PetscErrorCode ierr;
661273d9f13SBarry Smith 
662273d9f13SBarry Smith   PetscFunctionBegin;
663273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
664273d9f13SBarry Smith 
665b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
666273d9f13SBarry Smith   A->data = (void*)b;
667273d9f13SBarry Smith 
66826283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
66926283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
670273d9f13SBarry Smith 
671273d9f13SBarry Smith   b->ctx                 = 0;
672ef66eb69SBarry Smith   b->vshift              = 0.0;
673ef66eb69SBarry Smith   b->vscale              = 1.0;
6740c0fd78eSBarry Smith   b->managescalingshifts = PETSC_TRUE;
675273d9f13SBarry Smith   A->assembled           = PETSC_TRUE;
676210f0121SBarry Smith   A->preallocated        = PETSC_FALSE;
6772205254eSKarl Rupp 
67817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
679273d9f13SBarry Smith   PetscFunctionReturn(0);
680273d9f13SBarry Smith }
681e51e0e81SBarry Smith 
6824b828684SBarry Smith /*@C
683052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
684ff756334SLois Curfman McInnes    private data storage format.
685e51e0e81SBarry Smith 
686c7fcc2eaSBarry Smith   Collective on MPI_Comm
687c7fcc2eaSBarry Smith 
688e51e0e81SBarry Smith    Input Parameters:
689c7fcc2eaSBarry Smith +  comm - MPI communicator
690c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
691c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
692c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
693c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
694c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
695e51e0e81SBarry Smith 
696ff756334SLois Curfman McInnes    Output Parameter:
69744cd7ae7SLois Curfman McInnes .  A - the matrix
698e51e0e81SBarry Smith 
699ff2fd236SBarry Smith    Level: advanced
700ff2fd236SBarry Smith 
701f39d1f56SLois Curfman McInnes   Usage:
7027b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
703f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
704c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
705f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
706f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
707f39d1f56SLois Curfman McInnes 
708ff756334SLois Curfman McInnes    Notes:
709ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
710ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
711ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
712e51e0e81SBarry Smith 
71395452b02SPatrick Sanan    Fortran Notes:
71495452b02SPatrick Sanan     To use this from Fortran with a ctx you must write an interface definition for this
715daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
716daf670e6SBarry Smith     in as the ctx argument.
717f38a8d56SBarry Smith 
718f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
719f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
720645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
721645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
722645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
723645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
724645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
725645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
726645985a0SLois Curfman McInnes    For example,
727f39d1f56SLois Curfman McInnes 
728f39d1f56SLois Curfman McInnes $
729f39d1f56SLois Curfman McInnes $     Vec x, y
7307b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
731645985a0SLois Curfman McInnes $     Mat A
732f39d1f56SLois Curfman McInnes $
733c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
734c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
735f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
736c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
737c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
738c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
739645985a0SLois Curfman McInnes $     MatMult(A,x,y);
740645985a0SLois Curfman McInnes $     MatDestroy(A);
741f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
742645985a0SLois Curfman McInnes $
743e51e0e81SBarry Smith 
7449b53d762SBarry Smith 
7459b53d762SBarry Smith    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), and MatScale() internally so these
7469b53d762SBarry Smith    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.
7479b53d762SBarry Smith 
7489b53d762SBarry Smith 
7490c0fd78eSBarry Smith     For rectangular matrices do all the scalings and shifts make sense?
7500c0fd78eSBarry Smith 
75195452b02SPatrick Sanan     Developers Notes:
75295452b02SPatrick Sanan     Regarding shifting and scaling. The general form is
7530c0fd78eSBarry Smith 
7540c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
7550c0fd78eSBarry Smith 
7560c0fd78eSBarry Smith       The order you apply the operations is important. For example if you have a dshift then
7570c0fd78eSBarry Smith       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
7580c0fd78eSBarry Smith       you get s*vscale*A + diag(shift)
7590c0fd78eSBarry Smith 
7600c0fd78eSBarry Smith           A is the user provided function.
7610c0fd78eSBarry Smith 
762*ad2f5c3fSBarry Smith    KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for
763*ad2f5c3fSBarry Smith    for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger
764*ad2f5c3fSBarry Smith    an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat);
765*ad2f5c3fSBarry Smith    each time the MATSHELL matrix has changed.
766*ad2f5c3fSBarry Smith 
767*ad2f5c3fSBarry Smith    Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided
768*ad2f5c3fSBarry Smith    with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale().
769*ad2f5c3fSBarry Smith 
7700b627109SLois Curfman McInnes .keywords: matrix, shell, create
7710b627109SLois Curfman McInnes 
7720c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
773e51e0e81SBarry Smith @*/
7747087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
775e51e0e81SBarry Smith {
776dfbe8321SBarry Smith   PetscErrorCode ierr;
777ed3cc1f0SBarry Smith 
7783a40ed3dSBarry Smith   PetscFunctionBegin;
779f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
780f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
781273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
782273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
7830fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
784273d9f13SBarry Smith   PetscFunctionReturn(0);
785c7fcc2eaSBarry Smith }
786c7fcc2eaSBarry Smith 
787c6866cfdSSatish Balay /*@
788273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
789c7fcc2eaSBarry Smith 
7903f9fe445SBarry Smith    Logically Collective on Mat
791c7fcc2eaSBarry Smith 
792273d9f13SBarry Smith     Input Parameters:
793273d9f13SBarry Smith +   mat - the shell matrix
794273d9f13SBarry Smith -   ctx - the context
795273d9f13SBarry Smith 
796273d9f13SBarry Smith    Level: advanced
797273d9f13SBarry Smith 
79895452b02SPatrick Sanan    Fortran Notes:
79995452b02SPatrick Sanan     To use this from Fortran you must write a Fortran interface definition for this
800daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
801273d9f13SBarry Smith 
802273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
8030bc0a719SBarry Smith @*/
8047087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
805273d9f13SBarry Smith {
806273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
807dfbe8321SBarry Smith   PetscErrorCode ierr;
808ace3abfcSBarry Smith   PetscBool      flg;
809273d9f13SBarry Smith 
810273d9f13SBarry Smith   PetscFunctionBegin;
8110700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
812251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
813273d9f13SBarry Smith   if (flg) {
814273d9f13SBarry Smith     shell->ctx = ctx;
815ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
8163a40ed3dSBarry Smith   PetscFunctionReturn(0);
817e51e0e81SBarry Smith }
818e51e0e81SBarry Smith 
8190c0fd78eSBarry Smith /*@
8200c0fd78eSBarry Smith     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
8210c0fd78eSBarry Smith           after MatCreateShell()
8220c0fd78eSBarry Smith 
8230c0fd78eSBarry Smith    Logically Collective on Mat
8240c0fd78eSBarry Smith 
8250c0fd78eSBarry Smith     Input Parameter:
8260c0fd78eSBarry Smith .   mat - the shell matrix
8270c0fd78eSBarry Smith 
8280c0fd78eSBarry Smith   Level: advanced
8290c0fd78eSBarry Smith 
8300c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
8310c0fd78eSBarry Smith @*/
8320c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A)
8330c0fd78eSBarry Smith {
8340c0fd78eSBarry Smith   PetscErrorCode ierr;
835976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
8360c0fd78eSBarry Smith   PetscBool      flg;
8370c0fd78eSBarry Smith 
8380c0fd78eSBarry Smith   PetscFunctionBegin;
8390c0fd78eSBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
8400c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr);
8410c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
842976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)A->data;
8430c0fd78eSBarry Smith   shell->managescalingshifts = PETSC_FALSE;
844976e8c5aSLisandro Dalcin   A->ops->diagonalset   = NULL;
845976e8c5aSLisandro Dalcin   A->ops->diagonalscale = NULL;
846976e8c5aSLisandro Dalcin   A->ops->scale         = NULL;
847976e8c5aSLisandro Dalcin   A->ops->shift         = NULL;
848976e8c5aSLisandro Dalcin   A->ops->axpy          = NULL;
8490c0fd78eSBarry Smith   PetscFunctionReturn(0);
8500c0fd78eSBarry Smith }
8510c0fd78eSBarry Smith 
852c16cb8f2SBarry Smith /*@C
853f3b1f45cSBarry Smith     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
854f3b1f45cSBarry Smith 
855f3b1f45cSBarry Smith    Logically Collective on Mat
856f3b1f45cSBarry Smith 
857f3b1f45cSBarry Smith     Input Parameters:
858f3b1f45cSBarry Smith +   mat - the shell matrix
859f3b1f45cSBarry Smith .   f - the function
860f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
861f3b1f45cSBarry Smith -   ctx - an optional context for the function
862f3b1f45cSBarry Smith 
863f3b1f45cSBarry Smith    Output Parameter:
864f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
865f3b1f45cSBarry Smith 
866f3b1f45cSBarry Smith    Options Database:
867f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
868f3b1f45cSBarry Smith 
869f3b1f45cSBarry Smith    Level: advanced
870f3b1f45cSBarry Smith 
87195452b02SPatrick Sanan    Fortran Notes:
87295452b02SPatrick Sanan     Not supported from Fortran
873f3b1f45cSBarry Smith 
874f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
875f3b1f45cSBarry Smith @*/
876f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
877f3b1f45cSBarry Smith {
878f3b1f45cSBarry Smith   PetscErrorCode ierr;
879f3b1f45cSBarry Smith   PetscInt       m,n;
880f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
881f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
88274e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
883f3b1f45cSBarry Smith 
884f3b1f45cSBarry Smith   PetscFunctionBegin;
885f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
886f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr);
887f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
888f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr);
889f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
890f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
891f3b1f45cSBarry Smith 
892f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
893f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mat,&Dmat);CHKERRQ(ierr);
894f3b1f45cSBarry Smith 
895f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
896f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
897f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
898f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
899f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
900f3b1f45cSBarry Smith     flag = PETSC_FALSE;
901f3b1f45cSBarry Smith     if (v) {
902fc7aafd1SBarry 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);
903f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
904f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
905f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
906f3b1f45cSBarry Smith     }
907f3b1f45cSBarry Smith   } else if (v) {
908fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
909f3b1f45cSBarry Smith   }
910f3b1f45cSBarry Smith   if (flg) *flg = flag;
911f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
912f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
913f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
914f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
915f3b1f45cSBarry Smith   PetscFunctionReturn(0);
916f3b1f45cSBarry Smith }
917f3b1f45cSBarry Smith 
918f3b1f45cSBarry Smith /*@C
919f3b1f45cSBarry Smith     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.
920f3b1f45cSBarry Smith 
921f3b1f45cSBarry Smith    Logically Collective on Mat
922f3b1f45cSBarry Smith 
923f3b1f45cSBarry Smith     Input Parameters:
924f3b1f45cSBarry Smith +   mat - the shell matrix
925f3b1f45cSBarry Smith .   f - the function
926f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
927f3b1f45cSBarry Smith -   ctx - an optional context for the function
928f3b1f45cSBarry Smith 
929f3b1f45cSBarry Smith    Output Parameter:
930f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
931f3b1f45cSBarry Smith 
932f3b1f45cSBarry Smith    Options Database:
933f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
934f3b1f45cSBarry Smith 
935f3b1f45cSBarry Smith    Level: advanced
936f3b1f45cSBarry Smith 
93795452b02SPatrick Sanan    Fortran Notes:
93895452b02SPatrick Sanan     Not supported from Fortran
939f3b1f45cSBarry Smith 
940f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
941f3b1f45cSBarry Smith @*/
942f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
943f3b1f45cSBarry Smith {
944f3b1f45cSBarry Smith   PetscErrorCode ierr;
945f3b1f45cSBarry Smith   Vec            x,y,z;
946f3b1f45cSBarry Smith   PetscInt       m,n,M,N;
947f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
948f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
94974e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
950f3b1f45cSBarry Smith 
951f3b1f45cSBarry Smith   PetscFunctionBegin;
952f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
953f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
954f3b1f45cSBarry Smith   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
955f3b1f45cSBarry Smith   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
956f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
957f3b1f45cSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
958f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
959f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
960f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
961f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
962f3b1f45cSBarry Smith   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
963f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr);
964f3b1f45cSBarry Smith 
965f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
966f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
967f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
968f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
969f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
970f3b1f45cSBarry Smith     flag = PETSC_FALSE;
971f3b1f45cSBarry Smith     if (v) {
972fc7aafd1SBarry 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);
973f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
974f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
975f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
976f3b1f45cSBarry Smith     }
977f3b1f45cSBarry Smith   } else if (v) {
978fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
979f3b1f45cSBarry Smith   }
980f3b1f45cSBarry Smith   if (flg) *flg = flag;
981f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
982f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
983f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
984f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
985f3b1f45cSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
986f3b1f45cSBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
987f3b1f45cSBarry Smith   ierr = VecDestroy(&z);CHKERRQ(ierr);
988f3b1f45cSBarry Smith   PetscFunctionReturn(0);
989f3b1f45cSBarry Smith }
990f3b1f45cSBarry Smith 
991f3b1f45cSBarry Smith /*@C
9920c0fd78eSBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
993e51e0e81SBarry Smith 
9943f9fe445SBarry Smith    Logically Collective on Mat
995fee21e36SBarry Smith 
996c7fcc2eaSBarry Smith     Input Parameters:
997c7fcc2eaSBarry Smith +   mat - the shell matrix
998c7fcc2eaSBarry Smith .   op - the name of the operation
999c7fcc2eaSBarry Smith -   f - the function that provides the operation.
1000c7fcc2eaSBarry Smith 
100115091d37SBarry Smith    Level: advanced
100215091d37SBarry Smith 
1003fae171e0SBarry Smith     Usage:
1004e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
1005f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
1006c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
10070b627109SLois Curfman McInnes 
1008a62d957aSLois Curfman McInnes     Notes:
1009e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
10101c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
1011a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
10121c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
1013a62d957aSLois Curfman McInnes 
1014a4d85fd7SBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
1015deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
1016deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
1017deebb3c3SLois Curfman McInnes     routines, e.g.,
1018a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1019a62d957aSLois Curfman McInnes 
10204aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
10214aa34b0aSBarry Smith     nonzero on failure.
10224aa34b0aSBarry Smith 
1023a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
1024a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
1025a62d957aSLois Curfman McInnes     set by MatCreateShell().
1026a62d957aSLois Curfman McInnes 
102795452b02SPatrick Sanan     Fortran Notes:
102895452b02SPatrick Sanan     For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
1029501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
1030501d9185SBarry Smith 
10310c0fd78eSBarry Smith     Use MatSetOperation() to set an operation for any matrix type
10320c0fd78eSBarry Smith 
1033a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
1034a62d957aSLois Curfman McInnes 
10350c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts()
1036e51e0e81SBarry Smith @*/
10377087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
1038e51e0e81SBarry Smith {
1039ace3abfcSBarry Smith   PetscBool      flg;
1040976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
1041976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
1042273d9f13SBarry Smith 
10433a40ed3dSBarry Smith   PetscFunctionBegin;
10440700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
10450c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
10460c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1047976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)mat->data;
1048976e8c5aSLisandro Dalcin 
10495edf6516SJed Brown   switch (op) {
10505edf6516SJed Brown   case MATOP_DESTROY:
105128f357e3SAlex Fikl     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
10525edf6516SJed Brown     break;
1053976e8c5aSLisandro Dalcin   case MATOP_VIEW:
1054976e8c5aSLisandro Dalcin     if (!mat->ops->viewnative) {
1055976e8c5aSLisandro Dalcin       mat->ops->viewnative = mat->ops->view;
1056976e8c5aSLisandro Dalcin     }
1057976e8c5aSLisandro Dalcin     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
1058976e8c5aSLisandro Dalcin     break;
1059976e8c5aSLisandro Dalcin   case MATOP_COPY:
1060976e8c5aSLisandro Dalcin     shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f;
1061976e8c5aSLisandro Dalcin     break;
10626464bf51SAlex Fikl   case MATOP_DIAGONAL_SET:
10630c0fd78eSBarry Smith   case MATOP_DIAGONAL_SCALE:
10640c0fd78eSBarry Smith   case MATOP_SHIFT:
10650c0fd78eSBarry Smith   case MATOP_SCALE:
10669f137db4SBarry Smith   case MATOP_AXPY:
10670c0fd78eSBarry Smith     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
10680c0fd78eSBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
10696464bf51SAlex Fikl     break;
10700c0fd78eSBarry Smith   case MATOP_GET_DIAGONAL:
1071976e8c5aSLisandro Dalcin     if (shell->managescalingshifts) {
1072976e8c5aSLisandro Dalcin       shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1073976e8c5aSLisandro Dalcin       mat->ops->getdiagonal   = MatGetDiagonal_Shell;
1074976e8c5aSLisandro Dalcin     } else {
1075976e8c5aSLisandro Dalcin       shell->ops->getdiagonal = NULL;
1076976e8c5aSLisandro Dalcin       mat->ops->getdiagonal   = (PetscErrorCode (*)(Mat,Vec))f;
107740a5a6c4SBarry Smith     }
10785edf6516SJed Brown     break;
10795edf6516SJed Brown   case MATOP_MULT:
10809f137db4SBarry Smith     if (shell->managescalingshifts) {
10819f137db4SBarry Smith       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10829f137db4SBarry Smith       mat->ops->mult   = MatMult_Shell;
1083976e8c5aSLisandro Dalcin     } else {
1084976e8c5aSLisandro Dalcin       shell->ops->mult = NULL;
1085976e8c5aSLisandro Dalcin       mat->ops->mult   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1086976e8c5aSLisandro Dalcin     }
10875edf6516SJed Brown     break;
10885edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
10899f137db4SBarry Smith     if (shell->managescalingshifts) {
10909f137db4SBarry Smith       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10915259c5a4SBarry Smith       mat->ops->multtranspose   = MatMultTranspose_Shell;
1092976e8c5aSLisandro Dalcin     } else {
1093976e8c5aSLisandro Dalcin       shell->ops->multtranspose = NULL;
1094976e8c5aSLisandro Dalcin       mat->ops->multtranspose   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1095976e8c5aSLisandro Dalcin     }
10965edf6516SJed Brown     break;
10975edf6516SJed Brown   default:
10985edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
10995ab264f3SAlp Dener     break;
1100a62d957aSLois Curfman McInnes   }
11013a40ed3dSBarry Smith   PetscFunctionReturn(0);
1102e51e0e81SBarry Smith }
1103f0479e8cSBarry Smith 
1104d4bb536fSBarry Smith /*@C
1105d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
1106d4bb536fSBarry Smith 
1107c7fcc2eaSBarry Smith     Not Collective
1108c7fcc2eaSBarry Smith 
1109d4bb536fSBarry Smith     Input Parameters:
1110c7fcc2eaSBarry Smith +   mat - the shell matrix
1111c7fcc2eaSBarry Smith -   op - the name of the operation
1112d4bb536fSBarry Smith 
1113d4bb536fSBarry Smith     Output Parameter:
1114d4bb536fSBarry Smith .   f - the function that provides the operation.
1115d4bb536fSBarry Smith 
111615091d37SBarry Smith     Level: advanced
111715091d37SBarry Smith 
1118d4bb536fSBarry Smith     Notes:
1119e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
1120d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
1121d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
1122d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
1123d4bb536fSBarry Smith 
1124d4bb536fSBarry Smith     All user-provided functions have the same calling
1125d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
1126d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
1127d4bb536fSBarry Smith     routines, e.g.,
1128d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1129d4bb536fSBarry Smith 
1130d4bb536fSBarry Smith     Within each user-defined routine, the user should call
1131d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
1132d4bb536fSBarry Smith     set by MatCreateShell().
1133d4bb536fSBarry Smith 
1134d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
1135d4bb536fSBarry Smith 
1136ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1137d4bb536fSBarry Smith @*/
11387087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
1139d4bb536fSBarry Smith {
1140ace3abfcSBarry Smith   PetscBool      flg;
1141976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
1142976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
1143273d9f13SBarry Smith 
11443a40ed3dSBarry Smith   PetscFunctionBegin;
11450700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1146976e8c5aSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1147976e8c5aSLisandro Dalcin   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1148976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)mat->data;
1149976e8c5aSLisandro Dalcin 
115028f357e3SAlex Fikl   switch (op) {
115128f357e3SAlex Fikl   case MATOP_DESTROY:
115228f357e3SAlex Fikl     *f = (void (*)(void))shell->ops->destroy;
115328f357e3SAlex Fikl     break;
115428f357e3SAlex Fikl   case MATOP_VIEW:
115528f357e3SAlex Fikl     *f = (void (*)(void))mat->ops->view;
115628f357e3SAlex Fikl     break;
1157976e8c5aSLisandro Dalcin   case MATOP_COPY:
1158976e8c5aSLisandro Dalcin     *f = (void (*)(void))shell->ops->copy;
1159976e8c5aSLisandro Dalcin     break;
1160976e8c5aSLisandro Dalcin   case MATOP_DIAGONAL_SET:
1161976e8c5aSLisandro Dalcin   case MATOP_DIAGONAL_SCALE:
1162976e8c5aSLisandro Dalcin   case MATOP_SHIFT:
1163976e8c5aSLisandro Dalcin   case MATOP_SCALE:
1164976e8c5aSLisandro Dalcin   case MATOP_AXPY:
1165976e8c5aSLisandro Dalcin     *f = (((void (**)(void))mat->ops)[op]);
1166976e8c5aSLisandro Dalcin     break;
1167976e8c5aSLisandro Dalcin   case MATOP_GET_DIAGONAL:
1168976e8c5aSLisandro Dalcin     if (shell->ops->getdiagonal)
1169976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->getdiagonal;
1170976e8c5aSLisandro Dalcin     else
1171976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1172976e8c5aSLisandro Dalcin     break;
1173976e8c5aSLisandro Dalcin   case MATOP_MULT:
1174976e8c5aSLisandro Dalcin     if (shell->ops->mult)
1175976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->mult;
1176976e8c5aSLisandro Dalcin     else
1177976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1178976e8c5aSLisandro Dalcin     break;
1179976e8c5aSLisandro Dalcin   case MATOP_MULT_TRANSPOSE:
1180976e8c5aSLisandro Dalcin     if (shell->ops->multtranspose)
1181976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->multtranspose;
1182976e8c5aSLisandro Dalcin     else
1183976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1184976e8c5aSLisandro Dalcin     break;
118528f357e3SAlex Fikl   default:
1186c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
1187d4bb536fSBarry Smith   }
11883a40ed3dSBarry Smith   PetscFunctionReturn(0);
1189d4bb536fSBarry Smith }
1190