xref: /petsc/src/mat/impls/shell/shell.c (revision 976e8c5a1adf031be8d0c92383b998f5a99293f1)
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 {
11*976e8c5aSLisandro Dalcin   /*  3 */ PetscErrorCode (*mult)(Mat,Vec,Vec);
12*976e8c5aSLisandro Dalcin   /*  5 */ PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
13*976e8c5aSLisandro Dalcin   /* 17 */ PetscErrorCode (*getdiagonal)(Mat,Vec);
14*976e8c5aSLisandro Dalcin   /* 43 */ PetscErrorCode (*copy)(Mat,Mat,MatStructure);
15*976e8c5aSLisandro 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 
141daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
142daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
143a62d957aSLois Curfman McInnes 
144a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context
145a62d957aSLois Curfman McInnes 
146ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
1470bc0a719SBarry Smith @*/
1488fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
149b4fd4287SBarry Smith {
150dfbe8321SBarry Smith   PetscErrorCode ierr;
151ace3abfcSBarry Smith   PetscBool      flg;
152273d9f13SBarry Smith 
1533a40ed3dSBarry Smith   PetscFunctionBegin;
1540700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1554482741eSBarry Smith   PetscValidPointer(ctx,2);
156251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
157940183f3SBarry Smith   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
158ce94432eSBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
1593a40ed3dSBarry Smith   PetscFunctionReturn(0);
160b4fd4287SBarry Smith }
161b4fd4287SBarry Smith 
162dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
163e51e0e81SBarry Smith {
164dfbe8321SBarry Smith   PetscErrorCode ierr;
165bf0cc555SLisandro Dalcin   Mat_Shell      *shell = (Mat_Shell*)mat->data;
166ed3cc1f0SBarry Smith 
1673a40ed3dSBarry Smith   PetscFunctionBegin;
16828f357e3SAlex Fikl   if (shell->ops->destroy) {
16928f357e3SAlex Fikl     ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr);
170bf0cc555SLisandro Dalcin   }
1710c0fd78eSBarry Smith   ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
1720c0fd78eSBarry Smith   ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
1730c0fd78eSBarry Smith   ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
1748fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
1758fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
1765edf6516SJed Brown   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
1775edf6516SJed Brown   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
1789f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
179bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
1803a40ed3dSBarry Smith   PetscFunctionReturn(0);
181e51e0e81SBarry Smith }
182e51e0e81SBarry Smith 
1837fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
1847fabda3fSAlex Fikl {
18528f357e3SAlex Fikl   Mat_Shell       *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
1867fabda3fSAlex Fikl   PetscErrorCode  ierr;
187cb8c8a77SBarry Smith   PetscBool       matflg;
1887fabda3fSAlex Fikl 
1897fabda3fSAlex Fikl   PetscFunctionBegin;
19028f357e3SAlex Fikl   ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr);
191cb8c8a77SBarry Smith   if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL");
1927fabda3fSAlex Fikl 
193cb8c8a77SBarry Smith   ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
194cb8c8a77SBarry Smith   ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
1957fabda3fSAlex Fikl 
196cb8c8a77SBarry Smith   if (shellA->ops->copy) {
19728f357e3SAlex Fikl     ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr);
198cb8c8a77SBarry Smith   }
1997fabda3fSAlex Fikl   shellB->vscale = shellA->vscale;
2007fabda3fSAlex Fikl   shellB->vshift = shellA->vshift;
2010c0fd78eSBarry Smith   if (shellA->dshift) {
2020c0fd78eSBarry Smith     if (!shellB->dshift) {
2030c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr);
2047fabda3fSAlex Fikl     }
2050c0fd78eSBarry Smith     ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr);
2067fabda3fSAlex Fikl   } else {
2079f137db4SBarry Smith     ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr);
2087fabda3fSAlex Fikl   }
2090c0fd78eSBarry Smith   if (shellA->left) {
2100c0fd78eSBarry Smith     if (!shellB->left) {
2110c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr);
2127fabda3fSAlex Fikl     }
2130c0fd78eSBarry Smith     ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr);
2147fabda3fSAlex Fikl   } else {
2159f137db4SBarry Smith     ierr = VecDestroy(&shellB->left);CHKERRQ(ierr);
2167fabda3fSAlex Fikl   }
2170c0fd78eSBarry Smith   if (shellA->right) {
2180c0fd78eSBarry Smith     if (!shellB->right) {
2190c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr);
2207fabda3fSAlex Fikl     }
2210c0fd78eSBarry Smith     ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr);
2227fabda3fSAlex Fikl   } else {
2239f137db4SBarry Smith     ierr = VecDestroy(&shellB->right);CHKERRQ(ierr);
2247fabda3fSAlex Fikl   }
22540e381c3SBarry Smith   ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr);
22640e381c3SBarry Smith   if (shellA->axpy) {
22740e381c3SBarry Smith     ierr                 = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr);
22840e381c3SBarry Smith     shellB->axpy        = shellA->axpy;
22940e381c3SBarry Smith     shellB->axpy_vscale = shellA->axpy_vscale;
23040e381c3SBarry Smith   }
2317fabda3fSAlex Fikl   PetscFunctionReturn(0);
2327fabda3fSAlex Fikl }
2337fabda3fSAlex Fikl 
234cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M)
235cb8c8a77SBarry Smith {
236cb8c8a77SBarry Smith   PetscErrorCode ierr;
237cb8c8a77SBarry Smith   void           *ctx;
238cb8c8a77SBarry Smith 
239cb8c8a77SBarry Smith   PetscFunctionBegin;
240cb8c8a77SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
241cb8c8a77SBarry Smith   ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr);
242cb8c8a77SBarry Smith   ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
243cb8c8a77SBarry Smith   PetscFunctionReturn(0);
244cb8c8a77SBarry Smith }
245cb8c8a77SBarry Smith 
246dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
247ef66eb69SBarry Smith {
248ef66eb69SBarry Smith   Mat_Shell        *shell = (Mat_Shell*)A->data;
249dfbe8321SBarry Smith   PetscErrorCode   ierr;
25025578ef6SJed Brown   Vec              xx;
251e3f487b0SBarry Smith   PetscObjectState instate,outstate;
252ef66eb69SBarry Smith 
253ef66eb69SBarry Smith   PetscFunctionBegin;
2548fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
255e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
256*976e8c5aSLisandro Dalcin   if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL");
25728f357e3SAlex Fikl   ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr);
258e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
259e3f487b0SBarry Smith   if (instate == outstate) {
260e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
261e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
262e3f487b0SBarry Smith   }
2638fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2648fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
2659f137db4SBarry Smith 
2669f137db4SBarry Smith   if (shell->axpy) {
2679f137db4SBarry Smith     if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);}
2689f137db4SBarry Smith     ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr);
2699f137db4SBarry Smith     ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
2709f137db4SBarry Smith   }
271ef66eb69SBarry Smith   PetscFunctionReturn(0);
272ef66eb69SBarry Smith }
273ef66eb69SBarry Smith 
2745edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
2755edf6516SJed Brown {
2765edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
2775edf6516SJed Brown   PetscErrorCode ierr;
2785edf6516SJed Brown 
2795edf6516SJed Brown   PetscFunctionBegin;
2805edf6516SJed Brown   if (y == z) {
2815edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
2825edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
283b8430d49SBarry Smith     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
2845edf6516SJed Brown   } else {
2855edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
2865edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
2875edf6516SJed Brown   }
2885edf6516SJed Brown   PetscFunctionReturn(0);
2895edf6516SJed Brown }
2905edf6516SJed Brown 
29118be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
29218be62a5SSatish Balay {
29318be62a5SSatish Balay   Mat_Shell        *shell = (Mat_Shell*)A->data;
29418be62a5SSatish Balay   PetscErrorCode   ierr;
29525578ef6SJed Brown   Vec              xx;
296e3f487b0SBarry Smith   PetscObjectState instate,outstate;
29718be62a5SSatish Balay 
29818be62a5SSatish Balay   PetscFunctionBegin;
2998fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
300e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
3010c0fd78eSBarry Smith   if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL");
30228f357e3SAlex Fikl   ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr);
303e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
304e3f487b0SBarry Smith   if (instate == outstate) {
305e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
306e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
307e3f487b0SBarry Smith   }
3088fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
3098fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
31018be62a5SSatish Balay   PetscFunctionReturn(0);
31118be62a5SSatish Balay }
31218be62a5SSatish Balay 
3135edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3145edf6516SJed Brown {
3155edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3165edf6516SJed Brown   PetscErrorCode ierr;
3175edf6516SJed Brown 
3185edf6516SJed Brown   PetscFunctionBegin;
3195edf6516SJed Brown   if (y == z) {
3205edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
3215edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
3225edf6516SJed Brown     ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr);
3235edf6516SJed Brown   } else {
3245edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
3255edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3265edf6516SJed Brown   }
3275edf6516SJed Brown   PetscFunctionReturn(0);
3285edf6516SJed Brown }
3295edf6516SJed Brown 
3300c0fd78eSBarry Smith /*
3310c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
3320c0fd78eSBarry Smith */
33318be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
33418be62a5SSatish Balay {
33518be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
33618be62a5SSatish Balay   PetscErrorCode ierr;
33718be62a5SSatish Balay 
33818be62a5SSatish Balay   PetscFunctionBegin;
3390c0fd78eSBarry Smith   if (shell->ops->getdiagonal) {
34028f357e3SAlex Fikl     ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr);
3410c0fd78eSBarry Smith   } else {
3420c0fd78eSBarry Smith     ierr = VecSet(v,0.0);CHKERRQ(ierr);
3430c0fd78eSBarry Smith   }
34418be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
3458fe8eb27SJed Brown   if (shell->dshift) {
3460c0fd78eSBarry Smith     ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr);
34718be62a5SSatish Balay   }
3480c0fd78eSBarry Smith   ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
3498fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
3508fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
35181c02519SBarry Smith   if (shell->axpy) {
35281c02519SBarry Smith     if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);}
35381c02519SBarry Smith     ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr);
35481c02519SBarry Smith     ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
35581c02519SBarry Smith   }
35618be62a5SSatish Balay   PetscFunctionReturn(0);
35718be62a5SSatish Balay }
35818be62a5SSatish Balay 
359f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
360ef66eb69SBarry Smith {
361ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3628fe8eb27SJed Brown   PetscErrorCode ierr;
363b24ad042SBarry Smith 
364ef66eb69SBarry Smith   PetscFunctionBegin;
3650c0fd78eSBarry Smith   if (shell->left || shell->right) {
3668fe8eb27SJed Brown     if (!shell->dshift) {
3670c0fd78eSBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
3680c0fd78eSBarry Smith       ierr = VecSet(shell->dshift,a);CHKERRQ(ierr);
3690c0fd78eSBarry Smith     } else {
3700c0fd78eSBarry Smith       if (shell->left)  {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
3710c0fd78eSBarry Smith       if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
3729f137db4SBarry Smith       ierr = VecShift(shell->dshift,a);CHKERRQ(ierr);
3730c0fd78eSBarry Smith     }
3748fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
3758fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
3768fe8eb27SJed Brown   } else shell->vshift += a;
377ef66eb69SBarry Smith   PetscFunctionReturn(0);
378ef66eb69SBarry Smith }
379ef66eb69SBarry Smith 
3806464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
3816464bf51SAlex Fikl {
3826464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
3836464bf51SAlex Fikl   PetscErrorCode ierr;
3846464bf51SAlex Fikl 
3856464bf51SAlex Fikl   PetscFunctionBegin;
3860c0fd78eSBarry Smith   if (ins == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported with INSERT_VALUES");
3870c0fd78eSBarry Smith   if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);}
3880c0fd78eSBarry Smith   if (shell->left || shell->right) {
38992fabd1fSBarry Smith     if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);}
3909f137db4SBarry Smith     if (shell->left && shell->right)  {
3919f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
3929f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
3939f137db4SBarry Smith     } else if (shell->left) {
3949f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
3959f137db4SBarry Smith     } else {
3969f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr);
3979f137db4SBarry Smith     }
3989f137db4SBarry Smith     if (!shell->dshift) {
3999f137db4SBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
4009f137db4SBarry Smith       ierr = VecCopy(shell->dshift,shell->right_work);CHKERRQ(ierr);
4019f137db4SBarry Smith     } else {
4029f137db4SBarry Smith       ierr = VecAXPY(shell->dshift,1.0,shell->right_work);CHKERRQ(ierr);
4039f137db4SBarry Smith     }
4040c0fd78eSBarry Smith   } else {
4050c0fd78eSBarry Smith     ierr = VecAXPY(shell->dshift,1.0,D);CHKERRQ(ierr);
4066464bf51SAlex Fikl   }
4076464bf51SAlex Fikl   PetscFunctionReturn(0);
4086464bf51SAlex Fikl }
4096464bf51SAlex Fikl 
410f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
411ef66eb69SBarry Smith {
412ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4138fe8eb27SJed Brown   PetscErrorCode ierr;
414b24ad042SBarry Smith 
415ef66eb69SBarry Smith   PetscFunctionBegin;
416f4df32b1SMatthew Knepley   shell->vscale *= a;
4170c0fd78eSBarry Smith   shell->vshift *= a;
4182205254eSKarl Rupp   if (shell->dshift) {
4192205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
4200c0fd78eSBarry Smith   }
42181c02519SBarry Smith   shell->axpy_vscale *= a;
4228fe8eb27SJed Brown   PetscFunctionReturn(0);
42318be62a5SSatish Balay }
4248fe8eb27SJed Brown 
4258fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
4268fe8eb27SJed Brown {
4278fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4288fe8eb27SJed Brown   PetscErrorCode ierr;
4298fe8eb27SJed Brown 
4308fe8eb27SJed Brown   PetscFunctionBegin;
43181c02519SBarry Smith   if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation");
4328fe8eb27SJed Brown   if (left) {
4330c0fd78eSBarry Smith     if (!shell->left) {
4340c0fd78eSBarry Smith       ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr);
4358fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
4360c0fd78eSBarry Smith     } else {
4370c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
43818be62a5SSatish Balay     }
439ef66eb69SBarry Smith   }
4408fe8eb27SJed Brown   if (right) {
4410c0fd78eSBarry Smith     if (!shell->right) {
4420c0fd78eSBarry Smith       ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr);
4438fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
4440c0fd78eSBarry Smith     } else {
4450c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
4468fe8eb27SJed Brown     }
4478fe8eb27SJed Brown   }
448ef66eb69SBarry Smith   PetscFunctionReturn(0);
449ef66eb69SBarry Smith }
450ef66eb69SBarry Smith 
451dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
452ef66eb69SBarry Smith {
453ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4540c0fd78eSBarry Smith   PetscErrorCode ierr;
455ef66eb69SBarry Smith 
456ef66eb69SBarry Smith   PetscFunctionBegin;
4578fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
458ef66eb69SBarry Smith     shell->vshift = 0.0;
459ef66eb69SBarry Smith     shell->vscale = 1.0;
4600c0fd78eSBarry Smith     ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
4610c0fd78eSBarry Smith     ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
4620c0fd78eSBarry Smith     ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
46381c02519SBarry Smith     ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
464ef66eb69SBarry Smith   }
465ef66eb69SBarry Smith   PetscFunctionReturn(0);
466ef66eb69SBarry Smith }
467ef66eb69SBarry Smith 
468cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
469b951964fSBarry Smith 
4703b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
4713b49f96aSBarry Smith {
4723b49f96aSBarry Smith   PetscFunctionBegin;
4733b49f96aSBarry Smith   *missing = PETSC_FALSE;
4743b49f96aSBarry Smith   PetscFunctionReturn(0);
4753b49f96aSBarry Smith }
4763b49f96aSBarry Smith 
4779f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
4789f137db4SBarry Smith {
4799f137db4SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4809f137db4SBarry Smith   PetscErrorCode ierr;
4819f137db4SBarry Smith 
4829f137db4SBarry Smith   PetscFunctionBegin;
4839f137db4SBarry Smith   ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr);
4849f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
4859f137db4SBarry Smith   shell->axpy        = X;
4869f137db4SBarry Smith   shell->axpy_vscale = a;
4879f137db4SBarry Smith   PetscFunctionReturn(0);
4889f137db4SBarry Smith }
4899f137db4SBarry Smith 
49009dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
49120563c6bSBarry Smith                                        0,
49220563c6bSBarry Smith                                        0,
4939f137db4SBarry Smith                                        0,
4940c0fd78eSBarry Smith                                 /* 4*/ MatMultAdd_Shell,
4959f137db4SBarry Smith                                        0,
4960c0fd78eSBarry Smith                                        MatMultTransposeAdd_Shell,
497b951964fSBarry Smith                                        0,
498b951964fSBarry Smith                                        0,
499b951964fSBarry Smith                                        0,
50097304618SKris Buschelman                                 /*10*/ 0,
501b951964fSBarry Smith                                        0,
502b951964fSBarry Smith                                        0,
503b951964fSBarry Smith                                        0,
504b951964fSBarry Smith                                        0,
50597304618SKris Buschelman                                 /*15*/ 0,
506b951964fSBarry Smith                                        0,
5070c0fd78eSBarry Smith                                        MatGetDiagonal_Shell,
5088fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
509b951964fSBarry Smith                                        0,
51097304618SKris Buschelman                                 /*20*/ 0,
511ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
512b951964fSBarry Smith                                        0,
513b951964fSBarry Smith                                        0,
514d519adbfSMatthew Knepley                                 /*24*/ 0,
515b951964fSBarry Smith                                        0,
516b951964fSBarry Smith                                        0,
517b951964fSBarry Smith                                        0,
518b951964fSBarry Smith                                        0,
519d519adbfSMatthew Knepley                                 /*29*/ 0,
520b951964fSBarry Smith                                        0,
521273d9f13SBarry Smith                                        0,
522b951964fSBarry Smith                                        0,
523b951964fSBarry Smith                                        0,
524cb8c8a77SBarry Smith                                 /*34*/ MatDuplicate_Shell,
525b951964fSBarry Smith                                        0,
526b951964fSBarry Smith                                        0,
52709dc0095SBarry Smith                                        0,
52809dc0095SBarry Smith                                        0,
5299f137db4SBarry Smith                                 /*39*/ MatAXPY_Shell,
53009dc0095SBarry Smith                                        0,
53109dc0095SBarry Smith                                        0,
53209dc0095SBarry Smith                                        0,
533cb8c8a77SBarry Smith                                        MatCopy_Shell,
534d519adbfSMatthew Knepley                                 /*44*/ 0,
535ef66eb69SBarry Smith                                        MatScale_Shell,
536ef66eb69SBarry Smith                                        MatShift_Shell,
5376464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
53809dc0095SBarry Smith                                        0,
539f73d5cc4SBarry Smith                                 /*49*/ 0,
54009dc0095SBarry Smith                                        0,
54109dc0095SBarry Smith                                        0,
54209dc0095SBarry Smith                                        0,
54309dc0095SBarry Smith                                        0,
544d519adbfSMatthew Knepley                                 /*54*/ 0,
54509dc0095SBarry Smith                                        0,
54609dc0095SBarry Smith                                        0,
54709dc0095SBarry Smith                                        0,
54809dc0095SBarry Smith                                        0,
549d519adbfSMatthew Knepley                                 /*59*/ 0,
550b9b97703SBarry Smith                                        MatDestroy_Shell,
55109dc0095SBarry Smith                                        0,
552357abbc8SBarry Smith                                        0,
553273d9f13SBarry Smith                                        0,
554d519adbfSMatthew Knepley                                 /*64*/ 0,
555273d9f13SBarry Smith                                        0,
556273d9f13SBarry Smith                                        0,
557273d9f13SBarry Smith                                        0,
558273d9f13SBarry Smith                                        0,
559d519adbfSMatthew Knepley                                 /*69*/ 0,
560273d9f13SBarry Smith                                        0,
561c87e5d42SMatthew Knepley                                        MatConvert_Shell,
562273d9f13SBarry Smith                                        0,
56397304618SKris Buschelman                                        0,
564d519adbfSMatthew Knepley                                 /*74*/ 0,
56597304618SKris Buschelman                                        0,
56697304618SKris Buschelman                                        0,
56797304618SKris Buschelman                                        0,
56897304618SKris Buschelman                                        0,
569d519adbfSMatthew Knepley                                 /*79*/ 0,
57097304618SKris Buschelman                                        0,
57197304618SKris Buschelman                                        0,
57297304618SKris Buschelman                                        0,
573865e5f61SKris Buschelman                                        0,
574d519adbfSMatthew Knepley                                 /*84*/ 0,
575865e5f61SKris Buschelman                                        0,
576865e5f61SKris Buschelman                                        0,
577865e5f61SKris Buschelman                                        0,
578865e5f61SKris Buschelman                                        0,
579d519adbfSMatthew Knepley                                 /*89*/ 0,
580865e5f61SKris Buschelman                                        0,
581865e5f61SKris Buschelman                                        0,
582865e5f61SKris Buschelman                                        0,
583865e5f61SKris Buschelman                                        0,
584d519adbfSMatthew Knepley                                 /*94*/ 0,
585865e5f61SKris Buschelman                                        0,
586865e5f61SKris Buschelman                                        0,
5873964eb88SJed Brown                                        0,
5883964eb88SJed Brown                                        0,
5893964eb88SJed Brown                                 /*99*/ 0,
5903964eb88SJed Brown                                        0,
5913964eb88SJed Brown                                        0,
5923964eb88SJed Brown                                        0,
5933964eb88SJed Brown                                        0,
5943964eb88SJed Brown                                /*104*/ 0,
5953964eb88SJed Brown                                        0,
5963964eb88SJed Brown                                        0,
5973964eb88SJed Brown                                        0,
5983964eb88SJed Brown                                        0,
5993964eb88SJed Brown                                /*109*/ 0,
6003964eb88SJed Brown                                        0,
6013964eb88SJed Brown                                        0,
6023964eb88SJed Brown                                        0,
6033b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
6043964eb88SJed Brown                                /*114*/ 0,
6053964eb88SJed Brown                                        0,
6063964eb88SJed Brown                                        0,
6073964eb88SJed Brown                                        0,
6083964eb88SJed Brown                                        0,
6093964eb88SJed Brown                                /*119*/ 0,
6103964eb88SJed Brown                                        0,
6113964eb88SJed Brown                                        0,
6123964eb88SJed Brown                                        0,
6133964eb88SJed Brown                                        0,
6143964eb88SJed Brown                                /*124*/ 0,
6153964eb88SJed Brown                                        0,
6163964eb88SJed Brown                                        0,
6173964eb88SJed Brown                                        0,
6183964eb88SJed Brown                                        0,
6193964eb88SJed Brown                                /*129*/ 0,
6203964eb88SJed Brown                                        0,
6213964eb88SJed Brown                                        0,
6223964eb88SJed Brown                                        0,
6233964eb88SJed Brown                                        0,
6243964eb88SJed Brown                                /*134*/ 0,
6253964eb88SJed Brown                                        0,
6263964eb88SJed Brown                                        0,
6273964eb88SJed Brown                                        0,
6283964eb88SJed Brown                                        0,
6293964eb88SJed Brown                                /*139*/ 0,
630f9426fe0SMark Adams                                        0,
6313964eb88SJed Brown                                        0
6323964eb88SJed Brown };
633273d9f13SBarry Smith 
6340bad9183SKris Buschelman /*MC
635fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
6360bad9183SKris Buschelman 
6370bad9183SKris Buschelman   Level: advanced
6380bad9183SKris Buschelman 
6390c0fd78eSBarry Smith .seealso: MatCreateShell()
6400bad9183SKris Buschelman M*/
6410bad9183SKris Buschelman 
6428cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
643273d9f13SBarry Smith {
644273d9f13SBarry Smith   Mat_Shell      *b;
645dfbe8321SBarry Smith   PetscErrorCode ierr;
646273d9f13SBarry Smith 
647273d9f13SBarry Smith   PetscFunctionBegin;
648273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
649273d9f13SBarry Smith 
650b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
651273d9f13SBarry Smith   A->data = (void*)b;
652273d9f13SBarry Smith 
65326283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
65426283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
655273d9f13SBarry Smith 
656273d9f13SBarry Smith   b->ctx                 = 0;
657ef66eb69SBarry Smith   b->vshift              = 0.0;
658ef66eb69SBarry Smith   b->vscale              = 1.0;
6590c0fd78eSBarry Smith   b->managescalingshifts = PETSC_TRUE;
660273d9f13SBarry Smith   A->assembled           = PETSC_TRUE;
661210f0121SBarry Smith   A->preallocated        = PETSC_FALSE;
6622205254eSKarl Rupp 
66317667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
664273d9f13SBarry Smith   PetscFunctionReturn(0);
665273d9f13SBarry Smith }
666e51e0e81SBarry Smith 
6674b828684SBarry Smith /*@C
668052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
669ff756334SLois Curfman McInnes    private data storage format.
670e51e0e81SBarry Smith 
671c7fcc2eaSBarry Smith   Collective on MPI_Comm
672c7fcc2eaSBarry Smith 
673e51e0e81SBarry Smith    Input Parameters:
674c7fcc2eaSBarry Smith +  comm - MPI communicator
675c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
676c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
677c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
678c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
679c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
680e51e0e81SBarry Smith 
681ff756334SLois Curfman McInnes    Output Parameter:
68244cd7ae7SLois Curfman McInnes .  A - the matrix
683e51e0e81SBarry Smith 
684ff2fd236SBarry Smith    Level: advanced
685ff2fd236SBarry Smith 
686f39d1f56SLois Curfman McInnes   Usage:
6877b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
688f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
689c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
690f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
691f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
692f39d1f56SLois Curfman McInnes 
693ff756334SLois Curfman McInnes    Notes:
694ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
695ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
696ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
697e51e0e81SBarry Smith 
698daf670e6SBarry Smith    Fortran Notes: 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 
7350c0fd78eSBarry Smith     Developers Notes: Regarding shifting and scaling. The general form is
7360c0fd78eSBarry Smith 
7370c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
7380c0fd78eSBarry Smith 
7390c0fd78eSBarry Smith       The order you apply the operations is important. For example if you have a dshift then
7400c0fd78eSBarry Smith       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
7410c0fd78eSBarry Smith       you get s*vscale*A + diag(shift)
7420c0fd78eSBarry Smith 
7430c0fd78eSBarry Smith           A is the user provided function.
7440c0fd78eSBarry Smith 
7450b627109SLois Curfman McInnes .keywords: matrix, shell, create
7460b627109SLois Curfman McInnes 
7470c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
748e51e0e81SBarry Smith @*/
7497087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
750e51e0e81SBarry Smith {
751dfbe8321SBarry Smith   PetscErrorCode ierr;
752ed3cc1f0SBarry Smith 
7533a40ed3dSBarry Smith   PetscFunctionBegin;
754f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
755f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
756273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
757273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
7580fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
759273d9f13SBarry Smith   PetscFunctionReturn(0);
760c7fcc2eaSBarry Smith }
761c7fcc2eaSBarry Smith 
762c6866cfdSSatish Balay /*@
763273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
764c7fcc2eaSBarry Smith 
7653f9fe445SBarry Smith    Logically Collective on Mat
766c7fcc2eaSBarry Smith 
767273d9f13SBarry Smith     Input Parameters:
768273d9f13SBarry Smith +   mat - the shell matrix
769273d9f13SBarry Smith -   ctx - the context
770273d9f13SBarry Smith 
771273d9f13SBarry Smith    Level: advanced
772273d9f13SBarry Smith 
773daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
774daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
775273d9f13SBarry Smith 
776273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
7770bc0a719SBarry Smith @*/
7787087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
779273d9f13SBarry Smith {
780273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
781dfbe8321SBarry Smith   PetscErrorCode ierr;
782ace3abfcSBarry Smith   PetscBool      flg;
783273d9f13SBarry Smith 
784273d9f13SBarry Smith   PetscFunctionBegin;
7850700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
786251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
787273d9f13SBarry Smith   if (flg) {
788273d9f13SBarry Smith     shell->ctx = ctx;
789ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
7903a40ed3dSBarry Smith   PetscFunctionReturn(0);
791e51e0e81SBarry Smith }
792e51e0e81SBarry Smith 
7930c0fd78eSBarry Smith /*@
7940c0fd78eSBarry Smith     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
7950c0fd78eSBarry Smith           after MatCreateShell()
7960c0fd78eSBarry Smith 
7970c0fd78eSBarry Smith    Logically Collective on Mat
7980c0fd78eSBarry Smith 
7990c0fd78eSBarry Smith     Input Parameter:
8000c0fd78eSBarry Smith .   mat - the shell matrix
8010c0fd78eSBarry Smith 
8020c0fd78eSBarry Smith   Level: advanced
8030c0fd78eSBarry Smith 
8040c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
8050c0fd78eSBarry Smith @*/
8060c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A)
8070c0fd78eSBarry Smith {
8080c0fd78eSBarry Smith   PetscErrorCode ierr;
809*976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
8100c0fd78eSBarry Smith   PetscBool      flg;
8110c0fd78eSBarry Smith 
8120c0fd78eSBarry Smith   PetscFunctionBegin;
8130c0fd78eSBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
8140c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr);
8150c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
816*976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)A->data;
8170c0fd78eSBarry Smith   shell->managescalingshifts = PETSC_FALSE;
818*976e8c5aSLisandro Dalcin   A->ops->diagonalset   = NULL;
819*976e8c5aSLisandro Dalcin   A->ops->diagonalscale = NULL;
820*976e8c5aSLisandro Dalcin   A->ops->scale         = NULL;
821*976e8c5aSLisandro Dalcin   A->ops->shift         = NULL;
822*976e8c5aSLisandro Dalcin   A->ops->axpy          = NULL;
8230c0fd78eSBarry Smith   PetscFunctionReturn(0);
8240c0fd78eSBarry Smith }
8250c0fd78eSBarry Smith 
826c16cb8f2SBarry Smith /*@C
827f3b1f45cSBarry Smith     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
828f3b1f45cSBarry Smith 
829f3b1f45cSBarry Smith    Logically Collective on Mat
830f3b1f45cSBarry Smith 
831f3b1f45cSBarry Smith     Input Parameters:
832f3b1f45cSBarry Smith +   mat - the shell matrix
833f3b1f45cSBarry Smith .   f - the function
834f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
835f3b1f45cSBarry Smith -   ctx - an optional context for the function
836f3b1f45cSBarry Smith 
837f3b1f45cSBarry Smith    Output Parameter:
838f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
839f3b1f45cSBarry Smith 
840f3b1f45cSBarry Smith    Options Database:
841f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
842f3b1f45cSBarry Smith 
843f3b1f45cSBarry Smith    Level: advanced
844f3b1f45cSBarry Smith 
845f3b1f45cSBarry Smith    Fortran Notes: Not supported from Fortran
846f3b1f45cSBarry Smith 
847f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
848f3b1f45cSBarry Smith @*/
849f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
850f3b1f45cSBarry Smith {
851f3b1f45cSBarry Smith   PetscErrorCode ierr;
852f3b1f45cSBarry Smith   PetscInt       m,n;
853f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
854f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
85574e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
856f3b1f45cSBarry Smith 
857f3b1f45cSBarry Smith   PetscFunctionBegin;
858f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
859f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr);
860f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
861f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr);
862f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
863f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
864f3b1f45cSBarry Smith 
865f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
866f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mat,&Dmat);CHKERRQ(ierr);
867f3b1f45cSBarry Smith 
868f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
869f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
870f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
871f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
872f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
873f3b1f45cSBarry Smith     flag = PETSC_FALSE;
874f3b1f45cSBarry Smith     if (v) {
875fc7aafd1SBarry 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);
876f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
877f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
878f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
879f3b1f45cSBarry Smith     }
880f3b1f45cSBarry Smith   } else if (v) {
881fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
882f3b1f45cSBarry Smith   }
883f3b1f45cSBarry Smith   if (flg) *flg = flag;
884f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
885f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
886f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
887f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
888f3b1f45cSBarry Smith   PetscFunctionReturn(0);
889f3b1f45cSBarry Smith }
890f3b1f45cSBarry Smith 
891f3b1f45cSBarry Smith /*@C
892f3b1f45cSBarry Smith     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.
893f3b1f45cSBarry Smith 
894f3b1f45cSBarry Smith    Logically Collective on Mat
895f3b1f45cSBarry Smith 
896f3b1f45cSBarry Smith     Input Parameters:
897f3b1f45cSBarry Smith +   mat - the shell matrix
898f3b1f45cSBarry Smith .   f - the function
899f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
900f3b1f45cSBarry Smith -   ctx - an optional context for the function
901f3b1f45cSBarry Smith 
902f3b1f45cSBarry Smith    Output Parameter:
903f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
904f3b1f45cSBarry Smith 
905f3b1f45cSBarry Smith    Options Database:
906f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
907f3b1f45cSBarry Smith 
908f3b1f45cSBarry Smith    Level: advanced
909f3b1f45cSBarry Smith 
910f3b1f45cSBarry Smith    Fortran Notes: Not supported from Fortran
911f3b1f45cSBarry Smith 
912f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
913f3b1f45cSBarry Smith @*/
914f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
915f3b1f45cSBarry Smith {
916f3b1f45cSBarry Smith   PetscErrorCode ierr;
917f3b1f45cSBarry Smith   Vec            x,y,z;
918f3b1f45cSBarry Smith   PetscInt       m,n,M,N;
919f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
920f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
92174e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
922f3b1f45cSBarry Smith 
923f3b1f45cSBarry Smith   PetscFunctionBegin;
924f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
925f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
926f3b1f45cSBarry Smith   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
927f3b1f45cSBarry Smith   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
928f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
929f3b1f45cSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
930f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
931f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
932f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
933f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
934f3b1f45cSBarry Smith   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
935f3b1f45cSBarry Smith   ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr);
936f3b1f45cSBarry Smith 
937f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
938f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
939f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
940f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
941f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
942f3b1f45cSBarry Smith     flag = PETSC_FALSE;
943f3b1f45cSBarry Smith     if (v) {
944fc7aafd1SBarry 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);
945f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
946f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
947f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
948f3b1f45cSBarry Smith     }
949f3b1f45cSBarry Smith   } else if (v) {
950fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
951f3b1f45cSBarry Smith   }
952f3b1f45cSBarry Smith   if (flg) *flg = flag;
953f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
954f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
955f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
956f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
957f3b1f45cSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
958f3b1f45cSBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
959f3b1f45cSBarry Smith   ierr = VecDestroy(&z);CHKERRQ(ierr);
960f3b1f45cSBarry Smith   PetscFunctionReturn(0);
961f3b1f45cSBarry Smith }
962f3b1f45cSBarry Smith 
963f3b1f45cSBarry Smith /*@C
9640c0fd78eSBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
965e51e0e81SBarry Smith 
9663f9fe445SBarry Smith    Logically Collective on Mat
967fee21e36SBarry Smith 
968c7fcc2eaSBarry Smith     Input Parameters:
969c7fcc2eaSBarry Smith +   mat - the shell matrix
970c7fcc2eaSBarry Smith .   op - the name of the operation
971c7fcc2eaSBarry Smith -   f - the function that provides the operation.
972c7fcc2eaSBarry Smith 
97315091d37SBarry Smith    Level: advanced
97415091d37SBarry Smith 
975fae171e0SBarry Smith     Usage:
976e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
977f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
978c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
9790b627109SLois Curfman McInnes 
980a62d957aSLois Curfman McInnes     Notes:
981e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
9821c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
983a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
9841c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
985a62d957aSLois Curfman McInnes 
986a4d85fd7SBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
987deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
988deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
989deebb3c3SLois Curfman McInnes     routines, e.g.,
990a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
991a62d957aSLois Curfman McInnes 
9924aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
9934aa34b0aSBarry Smith     nonzero on failure.
9944aa34b0aSBarry Smith 
995a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
996a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
997a62d957aSLois Curfman McInnes     set by MatCreateShell().
998a62d957aSLois Curfman McInnes 
9992a7a6963SBarry Smith     Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
1000501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
1001501d9185SBarry Smith 
10020c0fd78eSBarry Smith     Use MatSetOperation() to set an operation for any matrix type
10030c0fd78eSBarry Smith 
1004a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
1005a62d957aSLois Curfman McInnes 
10060c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts()
1007e51e0e81SBarry Smith @*/
10087087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
1009e51e0e81SBarry Smith {
1010ace3abfcSBarry Smith   PetscBool      flg;
1011*976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
1012*976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
1013273d9f13SBarry Smith 
10143a40ed3dSBarry Smith   PetscFunctionBegin;
10150700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
10160c0fd78eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
10170c0fd78eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1018*976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)mat->data;
1019*976e8c5aSLisandro Dalcin 
10205edf6516SJed Brown   switch (op) {
10215edf6516SJed Brown   case MATOP_DESTROY:
102228f357e3SAlex Fikl     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
10235edf6516SJed Brown     break;
1024*976e8c5aSLisandro Dalcin   case MATOP_VIEW:
1025*976e8c5aSLisandro Dalcin     if (!mat->ops->viewnative) {
1026*976e8c5aSLisandro Dalcin       mat->ops->viewnative = mat->ops->view;
1027*976e8c5aSLisandro Dalcin     }
1028*976e8c5aSLisandro Dalcin     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
1029*976e8c5aSLisandro Dalcin     break;
1030*976e8c5aSLisandro Dalcin   case MATOP_COPY:
1031*976e8c5aSLisandro Dalcin     shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f;
1032*976e8c5aSLisandro Dalcin     break;
10336464bf51SAlex Fikl   case MATOP_DIAGONAL_SET:
10340c0fd78eSBarry Smith   case MATOP_DIAGONAL_SCALE:
10350c0fd78eSBarry Smith   case MATOP_SHIFT:
10360c0fd78eSBarry Smith   case MATOP_SCALE:
10379f137db4SBarry Smith   case MATOP_AXPY:
10380c0fd78eSBarry Smith     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
10390c0fd78eSBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
10406464bf51SAlex Fikl     break;
10410c0fd78eSBarry Smith   case MATOP_GET_DIAGONAL:
1042*976e8c5aSLisandro Dalcin     if (shell->managescalingshifts) {
1043*976e8c5aSLisandro Dalcin       shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1044*976e8c5aSLisandro Dalcin       mat->ops->getdiagonal   = MatGetDiagonal_Shell;
1045*976e8c5aSLisandro Dalcin     } else {
1046*976e8c5aSLisandro Dalcin       shell->ops->getdiagonal = NULL;
1047*976e8c5aSLisandro Dalcin       mat->ops->getdiagonal   = (PetscErrorCode (*)(Mat,Vec))f;
104840a5a6c4SBarry Smith     }
10495edf6516SJed Brown     break;
10505edf6516SJed Brown   case MATOP_MULT:
10519f137db4SBarry Smith     if (shell->managescalingshifts) {
10529f137db4SBarry Smith       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10539f137db4SBarry Smith       mat->ops->mult   = MatMult_Shell;
1054*976e8c5aSLisandro Dalcin     } else {
1055*976e8c5aSLisandro Dalcin       shell->ops->mult = NULL;
1056*976e8c5aSLisandro Dalcin       mat->ops->mult   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1057*976e8c5aSLisandro Dalcin     }
10585edf6516SJed Brown     break;
10595edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
10609f137db4SBarry Smith     if (shell->managescalingshifts) {
10619f137db4SBarry Smith       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
10625259c5a4SBarry Smith       mat->ops->multtranspose   = MatMultTranspose_Shell;
1063*976e8c5aSLisandro Dalcin     } else {
1064*976e8c5aSLisandro Dalcin       shell->ops->multtranspose = NULL;
1065*976e8c5aSLisandro Dalcin       mat->ops->multtranspose   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1066*976e8c5aSLisandro Dalcin     }
10675edf6516SJed Brown     break;
10685edf6516SJed Brown   default:
10695edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
10705ab264f3SAlp Dener     break;
1071a62d957aSLois Curfman McInnes   }
10723a40ed3dSBarry Smith   PetscFunctionReturn(0);
1073e51e0e81SBarry Smith }
1074f0479e8cSBarry Smith 
1075d4bb536fSBarry Smith /*@C
1076d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
1077d4bb536fSBarry Smith 
1078c7fcc2eaSBarry Smith     Not Collective
1079c7fcc2eaSBarry Smith 
1080d4bb536fSBarry Smith     Input Parameters:
1081c7fcc2eaSBarry Smith +   mat - the shell matrix
1082c7fcc2eaSBarry Smith -   op - the name of the operation
1083d4bb536fSBarry Smith 
1084d4bb536fSBarry Smith     Output Parameter:
1085d4bb536fSBarry Smith .   f - the function that provides the operation.
1086d4bb536fSBarry Smith 
108715091d37SBarry Smith     Level: advanced
108815091d37SBarry Smith 
1089d4bb536fSBarry Smith     Notes:
1090e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
1091d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
1092d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
1093d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
1094d4bb536fSBarry Smith 
1095d4bb536fSBarry Smith     All user-provided functions have the same calling
1096d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
1097d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
1098d4bb536fSBarry Smith     routines, e.g.,
1099d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1100d4bb536fSBarry Smith 
1101d4bb536fSBarry Smith     Within each user-defined routine, the user should call
1102d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
1103d4bb536fSBarry Smith     set by MatCreateShell().
1104d4bb536fSBarry Smith 
1105d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
1106d4bb536fSBarry Smith 
1107ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1108d4bb536fSBarry Smith @*/
11097087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
1110d4bb536fSBarry Smith {
1111ace3abfcSBarry Smith   PetscBool      flg;
1112*976e8c5aSLisandro Dalcin   Mat_Shell      *shell;
1113*976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
1114273d9f13SBarry Smith 
11153a40ed3dSBarry Smith   PetscFunctionBegin;
11160700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1117*976e8c5aSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1118*976e8c5aSLisandro Dalcin   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1119*976e8c5aSLisandro Dalcin   shell = (Mat_Shell*)mat->data;
1120*976e8c5aSLisandro Dalcin 
112128f357e3SAlex Fikl   switch (op) {
112228f357e3SAlex Fikl   case MATOP_DESTROY:
112328f357e3SAlex Fikl     *f = (void (*)(void))shell->ops->destroy;
112428f357e3SAlex Fikl     break;
112528f357e3SAlex Fikl   case MATOP_VIEW:
112628f357e3SAlex Fikl     *f = (void (*)(void))mat->ops->view;
112728f357e3SAlex Fikl     break;
1128*976e8c5aSLisandro Dalcin   case MATOP_COPY:
1129*976e8c5aSLisandro Dalcin     *f = (void (*)(void))shell->ops->copy;
1130*976e8c5aSLisandro Dalcin     break;
1131*976e8c5aSLisandro Dalcin   case MATOP_DIAGONAL_SET:
1132*976e8c5aSLisandro Dalcin   case MATOP_DIAGONAL_SCALE:
1133*976e8c5aSLisandro Dalcin   case MATOP_SHIFT:
1134*976e8c5aSLisandro Dalcin   case MATOP_SCALE:
1135*976e8c5aSLisandro Dalcin   case MATOP_AXPY:
1136*976e8c5aSLisandro Dalcin     *f = (((void (**)(void))mat->ops)[op]);
1137*976e8c5aSLisandro Dalcin     break;
1138*976e8c5aSLisandro Dalcin   case MATOP_GET_DIAGONAL:
1139*976e8c5aSLisandro Dalcin     if (shell->ops->getdiagonal)
1140*976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->getdiagonal;
1141*976e8c5aSLisandro Dalcin     else
1142*976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1143*976e8c5aSLisandro Dalcin     break;
1144*976e8c5aSLisandro Dalcin   case MATOP_MULT:
1145*976e8c5aSLisandro Dalcin     if (shell->ops->mult)
1146*976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->mult;
1147*976e8c5aSLisandro Dalcin     else
1148*976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1149*976e8c5aSLisandro Dalcin     break;
1150*976e8c5aSLisandro Dalcin   case MATOP_MULT_TRANSPOSE:
1151*976e8c5aSLisandro Dalcin     if (shell->ops->multtranspose)
1152*976e8c5aSLisandro Dalcin       *f = (void (*)(void))shell->ops->multtranspose;
1153*976e8c5aSLisandro Dalcin     else
1154*976e8c5aSLisandro Dalcin       *f = (((void (**)(void))mat->ops)[op]);
1155*976e8c5aSLisandro Dalcin     break;
115628f357e3SAlex Fikl   default:
1157c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
1158d4bb536fSBarry Smith   }
11593a40ed3dSBarry Smith   PetscFunctionReturn(0);
1160d4bb536fSBarry Smith }
1161