xref: /petsc/src/mat/impls/shell/shell.c (revision 00849d43f4426ca23107734077286c9795b00083)
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);
3235edf6516SJed Brown     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 
3796464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
3806464bf51SAlex Fikl {
3816464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
3826464bf51SAlex Fikl   PetscErrorCode ierr;
3836464bf51SAlex Fikl 
3846464bf51SAlex Fikl   PetscFunctionBegin;
3850c0fd78eSBarry Smith   if (ins == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported with INSERT_VALUES");
3860c0fd78eSBarry Smith   if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);}
3870c0fd78eSBarry Smith   if (shell->left || shell->right) {
38892fabd1fSBarry Smith     if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);}
3899f137db4SBarry Smith     if (shell->left && shell->right)  {
3909f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
3919f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
3929f137db4SBarry Smith     } else if (shell->left) {
3939f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
3949f137db4SBarry Smith     } else {
3959f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr);
3969f137db4SBarry Smith     }
3979f137db4SBarry Smith     if (!shell->dshift) {
3989f137db4SBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
3999f137db4SBarry Smith       ierr = VecCopy(shell->dshift,shell->right_work);CHKERRQ(ierr);
4009f137db4SBarry Smith     } else {
4019f137db4SBarry Smith       ierr = VecAXPY(shell->dshift,1.0,shell->right_work);CHKERRQ(ierr);
4029f137db4SBarry Smith     }
4030c0fd78eSBarry Smith   } else {
4040c0fd78eSBarry Smith     ierr = VecAXPY(shell->dshift,1.0,D);CHKERRQ(ierr);
4056464bf51SAlex Fikl   }
4066464bf51SAlex Fikl   PetscFunctionReturn(0);
4076464bf51SAlex Fikl }
4086464bf51SAlex Fikl 
409f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
410ef66eb69SBarry Smith {
411ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4128fe8eb27SJed Brown   PetscErrorCode ierr;
413b24ad042SBarry Smith 
414ef66eb69SBarry Smith   PetscFunctionBegin;
415f4df32b1SMatthew Knepley   shell->vscale *= a;
4160c0fd78eSBarry Smith   shell->vshift *= a;
4172205254eSKarl Rupp   if (shell->dshift) {
4182205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
4190c0fd78eSBarry Smith   }
42081c02519SBarry Smith   shell->axpy_vscale *= a;
4218fe8eb27SJed Brown   PetscFunctionReturn(0);
42218be62a5SSatish Balay }
4238fe8eb27SJed Brown 
4248fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
4258fe8eb27SJed Brown {
4268fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4278fe8eb27SJed Brown   PetscErrorCode ierr;
4288fe8eb27SJed Brown 
4298fe8eb27SJed Brown   PetscFunctionBegin;
43081c02519SBarry Smith   if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation");
4318fe8eb27SJed Brown   if (left) {
4320c0fd78eSBarry Smith     if (!shell->left) {
4330c0fd78eSBarry Smith       ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr);
4348fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
4350c0fd78eSBarry Smith     } else {
4360c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
43718be62a5SSatish Balay     }
438ef66eb69SBarry Smith   }
4398fe8eb27SJed Brown   if (right) {
4400c0fd78eSBarry Smith     if (!shell->right) {
4410c0fd78eSBarry Smith       ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr);
4428fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
4430c0fd78eSBarry Smith     } else {
4440c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
4458fe8eb27SJed Brown     }
4468fe8eb27SJed Brown   }
447ef66eb69SBarry Smith   PetscFunctionReturn(0);
448ef66eb69SBarry Smith }
449ef66eb69SBarry Smith 
450dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
451ef66eb69SBarry Smith {
452ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4530c0fd78eSBarry Smith   PetscErrorCode ierr;
454ef66eb69SBarry Smith 
455ef66eb69SBarry Smith   PetscFunctionBegin;
4568fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
457ef66eb69SBarry Smith     shell->vshift = 0.0;
458ef66eb69SBarry Smith     shell->vscale = 1.0;
4590c0fd78eSBarry Smith     ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
4600c0fd78eSBarry Smith     ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
4610c0fd78eSBarry Smith     ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
46281c02519SBarry Smith     ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
463ef66eb69SBarry Smith   }
464ef66eb69SBarry Smith   PetscFunctionReturn(0);
465ef66eb69SBarry Smith }
466ef66eb69SBarry Smith 
467cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
468b951964fSBarry Smith 
4693b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
4703b49f96aSBarry Smith {
4713b49f96aSBarry Smith   PetscFunctionBegin;
4723b49f96aSBarry Smith   *missing = PETSC_FALSE;
4733b49f96aSBarry Smith   PetscFunctionReturn(0);
4743b49f96aSBarry Smith }
4753b49f96aSBarry Smith 
4769f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
4779f137db4SBarry Smith {
4789f137db4SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4799f137db4SBarry Smith   PetscErrorCode ierr;
4809f137db4SBarry Smith 
4819f137db4SBarry Smith   PetscFunctionBegin;
4829f137db4SBarry Smith   ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr);
4839f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
4849f137db4SBarry Smith   shell->axpy        = X;
4859f137db4SBarry Smith   shell->axpy_vscale = a;
4869f137db4SBarry Smith   PetscFunctionReturn(0);
4879f137db4SBarry Smith }
4889f137db4SBarry Smith 
48909dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
49020563c6bSBarry Smith                                        0,
49120563c6bSBarry Smith                                        0,
4929f137db4SBarry Smith                                        0,
4930c0fd78eSBarry Smith                                 /* 4*/ MatMultAdd_Shell,
4949f137db4SBarry Smith                                        0,
4950c0fd78eSBarry Smith                                        MatMultTransposeAdd_Shell,
496b951964fSBarry Smith                                        0,
497b951964fSBarry Smith                                        0,
498b951964fSBarry Smith                                        0,
49997304618SKris Buschelman                                 /*10*/ 0,
500b951964fSBarry Smith                                        0,
501b951964fSBarry Smith                                        0,
502b951964fSBarry Smith                                        0,
503b951964fSBarry Smith                                        0,
50497304618SKris Buschelman                                 /*15*/ 0,
505b951964fSBarry Smith                                        0,
506*00849d43SBarry Smith                                        0,
5078fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
508b951964fSBarry Smith                                        0,
50997304618SKris Buschelman                                 /*20*/ 0,
510ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
511b951964fSBarry Smith                                        0,
512b951964fSBarry Smith                                        0,
513d519adbfSMatthew Knepley                                 /*24*/ 0,
514b951964fSBarry Smith                                        0,
515b951964fSBarry Smith                                        0,
516b951964fSBarry Smith                                        0,
517b951964fSBarry Smith                                        0,
518d519adbfSMatthew Knepley                                 /*29*/ 0,
519b951964fSBarry Smith                                        0,
520273d9f13SBarry Smith                                        0,
521b951964fSBarry Smith                                        0,
522b951964fSBarry Smith                                        0,
523cb8c8a77SBarry Smith                                 /*34*/ MatDuplicate_Shell,
524b951964fSBarry Smith                                        0,
525b951964fSBarry Smith                                        0,
52609dc0095SBarry Smith                                        0,
52709dc0095SBarry Smith                                        0,
5289f137db4SBarry Smith                                 /*39*/ MatAXPY_Shell,
52909dc0095SBarry Smith                                        0,
53009dc0095SBarry Smith                                        0,
53109dc0095SBarry Smith                                        0,
532cb8c8a77SBarry Smith                                        MatCopy_Shell,
533d519adbfSMatthew Knepley                                 /*44*/ 0,
534ef66eb69SBarry Smith                                        MatScale_Shell,
535ef66eb69SBarry Smith                                        MatShift_Shell,
5366464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
53709dc0095SBarry Smith                                        0,
538f73d5cc4SBarry Smith                                 /*49*/ 0,
53909dc0095SBarry Smith                                        0,
54009dc0095SBarry Smith                                        0,
54109dc0095SBarry Smith                                        0,
54209dc0095SBarry Smith                                        0,
543d519adbfSMatthew Knepley                                 /*54*/ 0,
54409dc0095SBarry Smith                                        0,
54509dc0095SBarry Smith                                        0,
54609dc0095SBarry Smith                                        0,
54709dc0095SBarry Smith                                        0,
548d519adbfSMatthew Knepley                                 /*59*/ 0,
549b9b97703SBarry Smith                                        MatDestroy_Shell,
55009dc0095SBarry Smith                                        0,
551357abbc8SBarry Smith                                        0,
552273d9f13SBarry Smith                                        0,
553d519adbfSMatthew Knepley                                 /*64*/ 0,
554273d9f13SBarry Smith                                        0,
555273d9f13SBarry Smith                                        0,
556273d9f13SBarry Smith                                        0,
557273d9f13SBarry Smith                                        0,
558d519adbfSMatthew Knepley                                 /*69*/ 0,
559273d9f13SBarry Smith                                        0,
560c87e5d42SMatthew Knepley                                        MatConvert_Shell,
561273d9f13SBarry Smith                                        0,
56297304618SKris Buschelman                                        0,
563d519adbfSMatthew Knepley                                 /*74*/ 0,
56497304618SKris Buschelman                                        0,
56597304618SKris Buschelman                                        0,
56697304618SKris Buschelman                                        0,
56797304618SKris Buschelman                                        0,
568d519adbfSMatthew Knepley                                 /*79*/ 0,
56997304618SKris Buschelman                                        0,
57097304618SKris Buschelman                                        0,
57197304618SKris Buschelman                                        0,
572865e5f61SKris Buschelman                                        0,
573d519adbfSMatthew Knepley                                 /*84*/ 0,
574865e5f61SKris Buschelman                                        0,
575865e5f61SKris Buschelman                                        0,
576865e5f61SKris Buschelman                                        0,
577865e5f61SKris Buschelman                                        0,
578d519adbfSMatthew Knepley                                 /*89*/ 0,
579865e5f61SKris Buschelman                                        0,
580865e5f61SKris Buschelman                                        0,
581865e5f61SKris Buschelman                                        0,
582865e5f61SKris Buschelman                                        0,
583d519adbfSMatthew Knepley                                 /*94*/ 0,
584865e5f61SKris Buschelman                                        0,
585865e5f61SKris Buschelman                                        0,
5863964eb88SJed Brown                                        0,
5873964eb88SJed Brown                                        0,
5883964eb88SJed Brown                                 /*99*/ 0,
5893964eb88SJed Brown                                        0,
5903964eb88SJed Brown                                        0,
5913964eb88SJed Brown                                        0,
5923964eb88SJed Brown                                        0,
5933964eb88SJed Brown                                /*104*/ 0,
5943964eb88SJed Brown                                        0,
5953964eb88SJed Brown                                        0,
5963964eb88SJed Brown                                        0,
5973964eb88SJed Brown                                        0,
5983964eb88SJed Brown                                /*109*/ 0,
5993964eb88SJed Brown                                        0,
6003964eb88SJed Brown                                        0,
6013964eb88SJed Brown                                        0,
6023b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
6033964eb88SJed Brown                                /*114*/ 0,
6043964eb88SJed Brown                                        0,
6053964eb88SJed Brown                                        0,
6063964eb88SJed Brown                                        0,
6073964eb88SJed Brown                                        0,
6083964eb88SJed Brown                                /*119*/ 0,
6093964eb88SJed Brown                                        0,
6103964eb88SJed Brown                                        0,
6113964eb88SJed Brown                                        0,
6123964eb88SJed Brown                                        0,
6133964eb88SJed Brown                                /*124*/ 0,
6143964eb88SJed Brown                                        0,
6153964eb88SJed Brown                                        0,
6163964eb88SJed Brown                                        0,
6173964eb88SJed Brown                                        0,
6183964eb88SJed Brown                                /*129*/ 0,
6193964eb88SJed Brown                                        0,
6203964eb88SJed Brown                                        0,
6213964eb88SJed Brown                                        0,
6223964eb88SJed Brown                                        0,
6233964eb88SJed Brown                                /*134*/ 0,
6243964eb88SJed Brown                                        0,
6253964eb88SJed Brown                                        0,
6263964eb88SJed Brown                                        0,
6273964eb88SJed Brown                                        0,
6283964eb88SJed Brown                                /*139*/ 0,
629f9426fe0SMark Adams                                        0,
6303964eb88SJed Brown                                        0
6313964eb88SJed Brown };
632273d9f13SBarry Smith 
6330bad9183SKris Buschelman /*MC
634fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
6350bad9183SKris Buschelman 
6360bad9183SKris Buschelman   Level: advanced
6370bad9183SKris Buschelman 
6380c0fd78eSBarry Smith .seealso: MatCreateShell()
6390bad9183SKris Buschelman M*/
6400bad9183SKris Buschelman 
6418cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
642273d9f13SBarry Smith {
643273d9f13SBarry Smith   Mat_Shell      *b;
644dfbe8321SBarry Smith   PetscErrorCode ierr;
645273d9f13SBarry Smith 
646273d9f13SBarry Smith   PetscFunctionBegin;
647273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
648273d9f13SBarry Smith 
649b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
650273d9f13SBarry Smith   A->data = (void*)b;
651273d9f13SBarry Smith 
65226283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
65326283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
654273d9f13SBarry Smith 
655273d9f13SBarry Smith   b->ctx                 = 0;
656ef66eb69SBarry Smith   b->vshift              = 0.0;
657ef66eb69SBarry Smith   b->vscale              = 1.0;
6580c0fd78eSBarry Smith   b->managescalingshifts = PETSC_TRUE;
659273d9f13SBarry Smith   A->assembled           = PETSC_TRUE;
660210f0121SBarry Smith   A->preallocated        = PETSC_FALSE;
6612205254eSKarl Rupp 
66217667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
663273d9f13SBarry Smith   PetscFunctionReturn(0);
664273d9f13SBarry Smith }
665e51e0e81SBarry Smith 
6664b828684SBarry Smith /*@C
667052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
668ff756334SLois Curfman McInnes    private data storage format.
669e51e0e81SBarry Smith 
670c7fcc2eaSBarry Smith   Collective on MPI_Comm
671c7fcc2eaSBarry Smith 
672e51e0e81SBarry Smith    Input Parameters:
673c7fcc2eaSBarry Smith +  comm - MPI communicator
674c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
675c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
676c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
677c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
678c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
679e51e0e81SBarry Smith 
680ff756334SLois Curfman McInnes    Output Parameter:
68144cd7ae7SLois Curfman McInnes .  A - the matrix
682e51e0e81SBarry Smith 
683ff2fd236SBarry Smith    Level: advanced
684ff2fd236SBarry Smith 
685f39d1f56SLois Curfman McInnes   Usage:
6867b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
687f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
688c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
689f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
690f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
691f39d1f56SLois Curfman McInnes 
692ff756334SLois Curfman McInnes    Notes:
693ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
694ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
695ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
696e51e0e81SBarry Smith 
69795452b02SPatrick Sanan    Fortran Notes:
69895452b02SPatrick Sanan     To use this from Fortran with a ctx you must write an interface definition for this
699daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
700daf670e6SBarry Smith     in as the ctx argument.
701f38a8d56SBarry Smith 
702f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
703f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
704645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
705645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
706645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
707645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
708645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
709645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
710645985a0SLois Curfman McInnes    For example,
711f39d1f56SLois Curfman McInnes 
712f39d1f56SLois Curfman McInnes $
713f39d1f56SLois Curfman McInnes $     Vec x, y
7147b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
715645985a0SLois Curfman McInnes $     Mat A
716f39d1f56SLois Curfman McInnes $
717c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
718c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
719f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
720c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
721c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
722c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
723645985a0SLois Curfman McInnes $     MatMult(A,x,y);
724645985a0SLois Curfman McInnes $     MatDestroy(A);
725f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
726645985a0SLois Curfman McInnes $
727e51e0e81SBarry Smith 
7289b53d762SBarry Smith 
7299b53d762SBarry Smith    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), and MatScale() internally so these
7309b53d762SBarry Smith    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.
7319b53d762SBarry Smith 
7329b53d762SBarry Smith 
7330c0fd78eSBarry Smith     For rectangular matrices do all the scalings and shifts make sense?
7340c0fd78eSBarry Smith 
73595452b02SPatrick Sanan     Developers Notes:
73695452b02SPatrick Sanan     Regarding shifting and scaling. The general form is
7370c0fd78eSBarry Smith 
7380c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
7390c0fd78eSBarry Smith 
7400c0fd78eSBarry Smith       The order you apply the operations is important. For example if you have a dshift then
7410c0fd78eSBarry Smith       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
7420c0fd78eSBarry Smith       you get s*vscale*A + diag(shift)
7430c0fd78eSBarry Smith 
7440c0fd78eSBarry Smith           A is the user provided function.
7450c0fd78eSBarry Smith 
7460b627109SLois Curfman McInnes .keywords: matrix, shell, create
7470b627109SLois Curfman McInnes 
7480c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
749e51e0e81SBarry Smith @*/
7507087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
751e51e0e81SBarry Smith {
752dfbe8321SBarry Smith   PetscErrorCode ierr;
753ed3cc1f0SBarry Smith 
7543a40ed3dSBarry Smith   PetscFunctionBegin;
755f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
756f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
757273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
758273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
7590fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
760273d9f13SBarry Smith   PetscFunctionReturn(0);
761c7fcc2eaSBarry Smith }
762c7fcc2eaSBarry Smith 
763c6866cfdSSatish Balay /*@
764273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
765c7fcc2eaSBarry Smith 
7663f9fe445SBarry Smith    Logically Collective on Mat
767c7fcc2eaSBarry Smith 
768273d9f13SBarry Smith     Input Parameters:
769273d9f13SBarry Smith +   mat - the shell matrix
770273d9f13SBarry Smith -   ctx - the context
771273d9f13SBarry Smith 
772273d9f13SBarry Smith    Level: advanced
773273d9f13SBarry Smith 
77495452b02SPatrick Sanan    Fortran Notes:
77595452b02SPatrick Sanan     To use this from Fortran you must write a Fortran interface definition for this
776daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
777273d9f13SBarry Smith 
778273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
7790bc0a719SBarry Smith @*/
7807087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
781273d9f13SBarry Smith {
782273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
783dfbe8321SBarry Smith   PetscErrorCode ierr;
784ace3abfcSBarry Smith   PetscBool      flg;
785273d9f13SBarry Smith 
786273d9f13SBarry Smith   PetscFunctionBegin;
7870700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
788251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
789273d9f13SBarry Smith   if (flg) {
790273d9f13SBarry Smith     shell->ctx = ctx;
791ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
7923a40ed3dSBarry Smith   PetscFunctionReturn(0);
793e51e0e81SBarry Smith }
794e51e0e81SBarry Smith 
7950c0fd78eSBarry Smith /*@
7960c0fd78eSBarry Smith     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
7970c0fd78eSBarry Smith           after MatCreateShell()
7980c0fd78eSBarry Smith 
7990c0fd78eSBarry Smith    Logically Collective on Mat
8000c0fd78eSBarry Smith 
8010c0fd78eSBarry Smith     Input Parameter:
8020c0fd78eSBarry Smith .   mat - the shell matrix
8030c0fd78eSBarry Smith 
8040c0fd78eSBarry Smith   Level: advanced
8050c0fd78eSBarry Smith 
8060c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
8070c0fd78eSBarry Smith @*/
8080c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A)
8090c0fd78eSBarry Smith {
8100c0fd78eSBarry Smith   PetscErrorCode ierr;
811976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
8120c0fd78eSBarry Smith   PetscBool      flg;
8130c0fd78eSBarry Smith 
8140c0fd78eSBarry Smith   PetscFunctionBegin;
8150c0fd78eSBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
8160c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr);
8170c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
818976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)A->data;
8190c0fd78eSBarry Smith   shell->managescalingshifts = PETSC_FALSE;
820976e8c5aSLisandro Dalcin   A->ops->diagonalset   = NULL;
821976e8c5aSLisandro Dalcin   A->ops->diagonalscale = NULL;
822976e8c5aSLisandro Dalcin   A->ops->scale         = NULL;
823976e8c5aSLisandro Dalcin   A->ops->shift         = NULL;
824976e8c5aSLisandro Dalcin   A->ops->axpy          = NULL;
8250c0fd78eSBarry Smith   PetscFunctionReturn(0);
8260c0fd78eSBarry Smith }
8270c0fd78eSBarry Smith 
828c16cb8f2SBarry Smith /*@C
829f3b1f45cSBarry Smith     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
830f3b1f45cSBarry Smith 
831f3b1f45cSBarry Smith    Logically Collective on Mat
832f3b1f45cSBarry Smith 
833f3b1f45cSBarry Smith     Input Parameters:
834f3b1f45cSBarry Smith +   mat - the shell matrix
835f3b1f45cSBarry Smith .   f - the function
836f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
837f3b1f45cSBarry Smith -   ctx - an optional context for the function
838f3b1f45cSBarry Smith 
839f3b1f45cSBarry Smith    Output Parameter:
840f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
841f3b1f45cSBarry Smith 
842f3b1f45cSBarry Smith    Options Database:
843f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
844f3b1f45cSBarry Smith 
845f3b1f45cSBarry Smith    Level: advanced
846f3b1f45cSBarry Smith 
84795452b02SPatrick Sanan    Fortran Notes:
84895452b02SPatrick Sanan     Not supported from Fortran
849f3b1f45cSBarry Smith 
850f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
851f3b1f45cSBarry Smith @*/
852f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
853f3b1f45cSBarry Smith {
854f3b1f45cSBarry Smith   PetscErrorCode ierr;
855f3b1f45cSBarry Smith   PetscInt       m,n;
856f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
857f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
85874e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
859f3b1f45cSBarry Smith 
860f3b1f45cSBarry Smith   PetscFunctionBegin;
861f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
862f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr);
863f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
864f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr);
865f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
866f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
867f3b1f45cSBarry Smith 
868f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
869f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mat,&Dmat);CHKERRQ(ierr);
870f3b1f45cSBarry Smith 
871f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
872f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
873f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
874f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
875f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
876f3b1f45cSBarry Smith     flag = PETSC_FALSE;
877f3b1f45cSBarry Smith     if (v) {
878fc7aafd1SBarry 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);
879f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
880f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
881f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
882f3b1f45cSBarry Smith     }
883f3b1f45cSBarry Smith   } else if (v) {
884fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
885f3b1f45cSBarry Smith   }
886f3b1f45cSBarry Smith   if (flg) *flg = flag;
887f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
888f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
889f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
890f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
891f3b1f45cSBarry Smith   PetscFunctionReturn(0);
892f3b1f45cSBarry Smith }
893f3b1f45cSBarry Smith 
894f3b1f45cSBarry Smith /*@C
895f3b1f45cSBarry Smith     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.
896f3b1f45cSBarry Smith 
897f3b1f45cSBarry Smith    Logically Collective on Mat
898f3b1f45cSBarry Smith 
899f3b1f45cSBarry Smith     Input Parameters:
900f3b1f45cSBarry Smith +   mat - the shell matrix
901f3b1f45cSBarry Smith .   f - the function
902f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
903f3b1f45cSBarry Smith -   ctx - an optional context for the function
904f3b1f45cSBarry Smith 
905f3b1f45cSBarry Smith    Output Parameter:
906f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
907f3b1f45cSBarry Smith 
908f3b1f45cSBarry Smith    Options Database:
909f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
910f3b1f45cSBarry Smith 
911f3b1f45cSBarry Smith    Level: advanced
912f3b1f45cSBarry Smith 
91395452b02SPatrick Sanan    Fortran Notes:
91495452b02SPatrick Sanan     Not supported from Fortran
915f3b1f45cSBarry Smith 
916f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
917f3b1f45cSBarry Smith @*/
918f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
919f3b1f45cSBarry Smith {
920f3b1f45cSBarry Smith   PetscErrorCode ierr;
921f3b1f45cSBarry Smith   Vec            x,y,z;
922f3b1f45cSBarry Smith   PetscInt       m,n,M,N;
923f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
924f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
92574e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
926f3b1f45cSBarry Smith 
927f3b1f45cSBarry Smith   PetscFunctionBegin;
928f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
929f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
930f3b1f45cSBarry Smith   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
931f3b1f45cSBarry Smith   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
932f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
933f3b1f45cSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
934f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
935f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
936f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
937f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
938f3b1f45cSBarry Smith   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
939f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr);
940f3b1f45cSBarry Smith 
941f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
942f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
943f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
944f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
945f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
946f3b1f45cSBarry Smith     flag = PETSC_FALSE;
947f3b1f45cSBarry Smith     if (v) {
948fc7aafd1SBarry 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);
949f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
950f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
951f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
952f3b1f45cSBarry Smith     }
953f3b1f45cSBarry Smith   } else if (v) {
954fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
955f3b1f45cSBarry Smith   }
956f3b1f45cSBarry Smith   if (flg) *flg = flag;
957f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
958f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
959f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
960f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
961f3b1f45cSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
962f3b1f45cSBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
963f3b1f45cSBarry Smith   ierr = VecDestroy(&z);CHKERRQ(ierr);
964f3b1f45cSBarry Smith   PetscFunctionReturn(0);
965f3b1f45cSBarry Smith }
966f3b1f45cSBarry Smith 
967f3b1f45cSBarry Smith /*@C
9680c0fd78eSBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
969e51e0e81SBarry Smith 
9703f9fe445SBarry Smith    Logically Collective on Mat
971fee21e36SBarry Smith 
972c7fcc2eaSBarry Smith     Input Parameters:
973c7fcc2eaSBarry Smith +   mat - the shell matrix
974c7fcc2eaSBarry Smith .   op - the name of the operation
975c7fcc2eaSBarry Smith -   f - the function that provides the operation.
976c7fcc2eaSBarry Smith 
97715091d37SBarry Smith    Level: advanced
97815091d37SBarry Smith 
979fae171e0SBarry Smith     Usage:
980e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
981f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
982c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
9830b627109SLois Curfman McInnes 
984a62d957aSLois Curfman McInnes     Notes:
985e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
9861c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
987a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
9881c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
989a62d957aSLois Curfman McInnes 
990a4d85fd7SBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
991deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
992deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
993deebb3c3SLois Curfman McInnes     routines, e.g.,
994a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
995a62d957aSLois Curfman McInnes 
9964aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
9974aa34b0aSBarry Smith     nonzero on failure.
9984aa34b0aSBarry Smith 
999a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
1000a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
1001a62d957aSLois Curfman McInnes     set by MatCreateShell().
1002a62d957aSLois Curfman McInnes 
100395452b02SPatrick Sanan     Fortran Notes:
100495452b02SPatrick Sanan     For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
1005501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
1006501d9185SBarry Smith 
10070c0fd78eSBarry Smith     Use MatSetOperation() to set an operation for any matrix type
10080c0fd78eSBarry Smith 
1009a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
1010a62d957aSLois Curfman McInnes 
10110c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts()
1012e51e0e81SBarry Smith @*/
10137087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
1014e51e0e81SBarry Smith {
1015ace3abfcSBarry Smith   PetscBool      flg;
1016976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
1017976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
1018273d9f13SBarry Smith 
10193a40ed3dSBarry Smith   PetscFunctionBegin;
10200700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
10210c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
10220c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1023976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)mat->data;
1024976e8c5aSLisandro Dalcin 
10255edf6516SJed Brown   switch (op) {
10265edf6516SJed Brown   case MATOP_DESTROY:
102728f357e3SAlex Fikl     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
10285edf6516SJed Brown     break;
1029976e8c5aSLisandro Dalcin   case MATOP_VIEW:
1030976e8c5aSLisandro Dalcin     if (!mat->ops->viewnative) {
1031976e8c5aSLisandro Dalcin       mat->ops->viewnative = mat->ops->view;
1032976e8c5aSLisandro Dalcin     }
1033976e8c5aSLisandro Dalcin     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
1034976e8c5aSLisandro Dalcin     break;
1035976e8c5aSLisandro Dalcin   case MATOP_COPY:
1036976e8c5aSLisandro Dalcin     shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f;
1037976e8c5aSLisandro Dalcin     break;
10386464bf51SAlex Fikl   case MATOP_DIAGONAL_SET:
10390c0fd78eSBarry Smith   case MATOP_DIAGONAL_SCALE:
10400c0fd78eSBarry Smith   case MATOP_SHIFT:
10410c0fd78eSBarry Smith   case MATOP_SCALE:
10429f137db4SBarry Smith   case MATOP_AXPY:
10430c0fd78eSBarry Smith     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
10440c0fd78eSBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
10456464bf51SAlex Fikl     break;
10460c0fd78eSBarry Smith   case MATOP_GET_DIAGONAL:
1047976e8c5aSLisandro Dalcin     if (shell->managescalingshifts) {
1048976e8c5aSLisandro Dalcin       shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1049976e8c5aSLisandro Dalcin       mat->ops->getdiagonal   = MatGetDiagonal_Shell;
1050976e8c5aSLisandro Dalcin     } else {
1051976e8c5aSLisandro Dalcin       shell->ops->getdiagonal = NULL;
1052976e8c5aSLisandro Dalcin       mat->ops->getdiagonal   = (PetscErrorCode (*)(Mat,Vec))f;
105340a5a6c4SBarry Smith     }
10545edf6516SJed Brown     break;
10555edf6516SJed Brown   case MATOP_MULT:
10569f137db4SBarry Smith     if (shell->managescalingshifts) {
10579f137db4SBarry Smith       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10589f137db4SBarry Smith       mat->ops->mult   = MatMult_Shell;
1059976e8c5aSLisandro Dalcin     } else {
1060976e8c5aSLisandro Dalcin       shell->ops->mult = NULL;
1061976e8c5aSLisandro Dalcin       mat->ops->mult   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1062976e8c5aSLisandro Dalcin     }
10635edf6516SJed Brown     break;
10645edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
10659f137db4SBarry Smith     if (shell->managescalingshifts) {
10669f137db4SBarry Smith       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10675259c5a4SBarry Smith       mat->ops->multtranspose   = MatMultTranspose_Shell;
1068976e8c5aSLisandro Dalcin     } else {
1069976e8c5aSLisandro Dalcin       shell->ops->multtranspose = NULL;
1070976e8c5aSLisandro Dalcin       mat->ops->multtranspose   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1071976e8c5aSLisandro Dalcin     }
10725edf6516SJed Brown     break;
10735edf6516SJed Brown   default:
10745edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
10755ab264f3SAlp Dener     break;
1076a62d957aSLois Curfman McInnes   }
10773a40ed3dSBarry Smith   PetscFunctionReturn(0);
1078e51e0e81SBarry Smith }
1079f0479e8cSBarry Smith 
1080d4bb536fSBarry Smith /*@C
1081d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
1082d4bb536fSBarry Smith 
1083c7fcc2eaSBarry Smith     Not Collective
1084c7fcc2eaSBarry Smith 
1085d4bb536fSBarry Smith     Input Parameters:
1086c7fcc2eaSBarry Smith +   mat - the shell matrix
1087c7fcc2eaSBarry Smith -   op - the name of the operation
1088d4bb536fSBarry Smith 
1089d4bb536fSBarry Smith     Output Parameter:
1090d4bb536fSBarry Smith .   f - the function that provides the operation.
1091d4bb536fSBarry Smith 
109215091d37SBarry Smith     Level: advanced
109315091d37SBarry Smith 
1094d4bb536fSBarry Smith     Notes:
1095e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
1096d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
1097d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
1098d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
1099d4bb536fSBarry Smith 
1100d4bb536fSBarry Smith     All user-provided functions have the same calling
1101d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
1102d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
1103d4bb536fSBarry Smith     routines, e.g.,
1104d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1105d4bb536fSBarry Smith 
1106d4bb536fSBarry Smith     Within each user-defined routine, the user should call
1107d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
1108d4bb536fSBarry Smith     set by MatCreateShell().
1109d4bb536fSBarry Smith 
1110d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
1111d4bb536fSBarry Smith 
1112ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1113d4bb536fSBarry Smith @*/
11147087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
1115d4bb536fSBarry Smith {
1116ace3abfcSBarry Smith   PetscBool      flg;
1117976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
1118976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
1119273d9f13SBarry Smith 
11203a40ed3dSBarry Smith   PetscFunctionBegin;
11210700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1122976e8c5aSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1123976e8c5aSLisandro Dalcin   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1124976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)mat->data;
1125976e8c5aSLisandro Dalcin 
112628f357e3SAlex Fikl   switch (op) {
112728f357e3SAlex Fikl   case MATOP_DESTROY:
112828f357e3SAlex Fikl     *f = (void (*)(void))shell->ops->destroy;
112928f357e3SAlex Fikl     break;
113028f357e3SAlex Fikl   case MATOP_VIEW:
113128f357e3SAlex Fikl     *f = (void (*)(void))mat->ops->view;
113228f357e3SAlex Fikl     break;
1133976e8c5aSLisandro Dalcin   case MATOP_COPY:
1134976e8c5aSLisandro Dalcin     *f = (void (*)(void))shell->ops->copy;
1135976e8c5aSLisandro Dalcin     break;
1136976e8c5aSLisandro Dalcin   case MATOP_DIAGONAL_SET:
1137976e8c5aSLisandro Dalcin   case MATOP_DIAGONAL_SCALE:
1138976e8c5aSLisandro Dalcin   case MATOP_SHIFT:
1139976e8c5aSLisandro Dalcin   case MATOP_SCALE:
1140976e8c5aSLisandro Dalcin   case MATOP_AXPY:
1141976e8c5aSLisandro Dalcin     *f = (((void (**)(void))mat->ops)[op]);
1142976e8c5aSLisandro Dalcin     break;
1143976e8c5aSLisandro Dalcin   case MATOP_GET_DIAGONAL:
1144976e8c5aSLisandro Dalcin     if (shell->ops->getdiagonal)
1145976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->getdiagonal;
1146976e8c5aSLisandro Dalcin     else
1147976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1148976e8c5aSLisandro Dalcin     break;
1149976e8c5aSLisandro Dalcin   case MATOP_MULT:
1150976e8c5aSLisandro Dalcin     if (shell->ops->mult)
1151976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->mult;
1152976e8c5aSLisandro Dalcin     else
1153976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1154976e8c5aSLisandro Dalcin     break;
1155976e8c5aSLisandro Dalcin   case MATOP_MULT_TRANSPOSE:
1156976e8c5aSLisandro Dalcin     if (shell->ops->multtranspose)
1157976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->multtranspose;
1158976e8c5aSLisandro Dalcin     else
1159976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1160976e8c5aSLisandro Dalcin     break;
116128f357e3SAlex Fikl   default:
1162c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
1163d4bb536fSBarry Smith   }
11643a40ed3dSBarry Smith   PetscFunctionReturn(0);
1165d4bb536fSBarry Smith }
1166