xref: /petsc/src/mat/impls/shell/shell.c (revision 5edf65166e22df6feda2ed6bd8a6e348419dee17)
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 
8b45d2f2cSJed Brown #include <petsc-private/matimpl.h>        /*I "petscmat.h" I*/
9b45d2f2cSJed Brown #include <petsc-private/vecimpl.h>
10e51e0e81SBarry Smith 
1120563c6bSBarry Smith typedef struct {
126849ba73SBarry Smith   PetscErrorCode (*destroy)(Mat);
136849ba73SBarry Smith   PetscErrorCode (*mult)(Mat,Vec,Vec);
1418be62a5SSatish Balay   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
1518be62a5SSatish Balay   PetscErrorCode (*getdiagonal)(Mat,Vec);
16ef66eb69SBarry Smith   PetscScalar    vscale,vshift;
178fe8eb27SJed Brown   Vec            dshift;
188fe8eb27SJed Brown   Vec            left,right;
198fe8eb27SJed Brown   Vec            dshift_owned,left_owned,right_owned;
208fe8eb27SJed Brown   Vec            left_work,right_work;
21*5edf6516SJed Brown   Vec            left_add_work,right_add_work;
228fe8eb27SJed Brown   PetscBool      usingscaled;
2320563c6bSBarry Smith   void           *ctx;
2488cf3e7dSBarry Smith } Mat_Shell;
258fe8eb27SJed Brown /*
268fe8eb27SJed Brown  The most general expression for the matrix is
278fe8eb27SJed Brown 
288fe8eb27SJed Brown  S = L (a A + B) R
298fe8eb27SJed Brown 
308fe8eb27SJed Brown  where
318fe8eb27SJed Brown  A is the matrix defined by the user's function
328fe8eb27SJed Brown  a is a scalar multiple
338fe8eb27SJed Brown  L is left scaling
348fe8eb27SJed Brown  R is right scaling
358fe8eb27SJed Brown  B is a diagonal shift defined by
368fe8eb27SJed Brown    diag(dshift) if the vector dshift is non-NULL
378fe8eb27SJed Brown    vscale*identity otherwise
388fe8eb27SJed Brown 
398fe8eb27SJed Brown  The following identities apply:
408fe8eb27SJed Brown 
418fe8eb27SJed Brown  Scale by c:
428fe8eb27SJed Brown   c [L (a A + B) R] = L [(a c) A + c B] R
438fe8eb27SJed Brown 
448fe8eb27SJed Brown  Shift by c:
458fe8eb27SJed Brown   [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R
468fe8eb27SJed Brown 
478fe8eb27SJed Brown  Diagonal scaling is achieved by simply multiplying with existing L and R vectors
488fe8eb27SJed Brown 
498fe8eb27SJed Brown  In the data structure:
508fe8eb27SJed Brown 
518fe8eb27SJed Brown  vscale=1.0  means no special scaling will be applied
528fe8eb27SJed Brown  dshift=NULL means a constant diagonal shift (fall back to vshift)
538fe8eb27SJed Brown  vshift=0.0  means no constant diagonal shift, note that vshift is only used if dshift is NULL
548fe8eb27SJed Brown */
558fe8eb27SJed Brown 
568fe8eb27SJed Brown static PetscErrorCode MatMult_Shell(Mat,Vec,Vec);
578fe8eb27SJed Brown static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec);
588fe8eb27SJed Brown static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec);
598fe8eb27SJed Brown 
608fe8eb27SJed Brown #undef __FUNCT__
618fe8eb27SJed Brown #define __FUNCT__ "MatShellUseScaledMethods"
628fe8eb27SJed Brown static PetscErrorCode MatShellUseScaledMethods(Mat Y)
638fe8eb27SJed Brown {
648fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)Y->data;
658fe8eb27SJed Brown 
668fe8eb27SJed Brown   PetscFunctionBegin;
678fe8eb27SJed Brown   if (shell->usingscaled) PetscFunctionReturn(0);
688fe8eb27SJed Brown   shell->mult  = Y->ops->mult;
698fe8eb27SJed Brown   Y->ops->mult = MatMult_Shell;
708fe8eb27SJed Brown   if (Y->ops->multtranspose) {
718fe8eb27SJed Brown     shell->multtranspose  = Y->ops->multtranspose;
728fe8eb27SJed Brown     Y->ops->multtranspose = MatMultTranspose_Shell;
738fe8eb27SJed Brown   }
748fe8eb27SJed Brown   if (Y->ops->getdiagonal) {
758fe8eb27SJed Brown     shell->getdiagonal  = Y->ops->getdiagonal;
768fe8eb27SJed Brown     Y->ops->getdiagonal = MatGetDiagonal_Shell;
778fe8eb27SJed Brown   }
788fe8eb27SJed Brown   shell->usingscaled = PETSC_TRUE;
798fe8eb27SJed Brown   PetscFunctionReturn(0);
808fe8eb27SJed Brown }
818fe8eb27SJed Brown 
828fe8eb27SJed Brown #undef __FUNCT__
838fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleLeft"
848fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
858fe8eb27SJed Brown {
868fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
878fe8eb27SJed Brown   PetscErrorCode ierr;
888fe8eb27SJed Brown 
898fe8eb27SJed Brown   PetscFunctionBegin;
90d63f877fSJed Brown   *xx = PETSC_NULL;
918fe8eb27SJed Brown   if (!shell->left) {
928fe8eb27SJed Brown     *xx = x;
938fe8eb27SJed Brown   } else {
948fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
958fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
968fe8eb27SJed Brown     *xx = shell->left_work;
978fe8eb27SJed Brown   }
988fe8eb27SJed Brown   PetscFunctionReturn(0);
998fe8eb27SJed Brown }
1008fe8eb27SJed Brown 
1018fe8eb27SJed Brown #undef __FUNCT__
1028fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleRight"
1038fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
1048fe8eb27SJed Brown {
1058fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1068fe8eb27SJed Brown   PetscErrorCode ierr;
1078fe8eb27SJed Brown 
1088fe8eb27SJed Brown   PetscFunctionBegin;
109d63f877fSJed Brown   *xx = PETSC_NULL;
1108fe8eb27SJed Brown   if (!shell->right) {
1118fe8eb27SJed Brown     *xx = x;
1128fe8eb27SJed Brown   } else {
1138fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
1148fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
1158fe8eb27SJed Brown     *xx = shell->right_work;
1168fe8eb27SJed Brown   }
1178fe8eb27SJed Brown   PetscFunctionReturn(0);
1188fe8eb27SJed Brown }
1198fe8eb27SJed Brown 
1208fe8eb27SJed Brown #undef __FUNCT__
1218fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleLeft"
1228fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
1238fe8eb27SJed Brown {
1248fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1258fe8eb27SJed Brown   PetscErrorCode ierr;
1268fe8eb27SJed Brown 
1278fe8eb27SJed Brown   PetscFunctionBegin;
1288fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
1298fe8eb27SJed Brown   PetscFunctionReturn(0);
1308fe8eb27SJed Brown }
1318fe8eb27SJed Brown 
1328fe8eb27SJed Brown #undef __FUNCT__
1338fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleRight"
1348fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
1358fe8eb27SJed Brown {
1368fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1378fe8eb27SJed Brown   PetscErrorCode ierr;
1388fe8eb27SJed Brown 
1398fe8eb27SJed Brown   PetscFunctionBegin;
1408fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
1418fe8eb27SJed Brown   PetscFunctionReturn(0);
1428fe8eb27SJed Brown }
1438fe8eb27SJed Brown 
1448fe8eb27SJed Brown #undef __FUNCT__
1458fe8eb27SJed Brown #define __FUNCT__ "MatShellShiftAndScale"
1468fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
1478fe8eb27SJed Brown {
1488fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1498fe8eb27SJed Brown   PetscErrorCode ierr;
1508fe8eb27SJed Brown 
1518fe8eb27SJed Brown   PetscFunctionBegin;
1528fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
1538fe8eb27SJed Brown     PetscInt          i,m;
1548fe8eb27SJed Brown     const PetscScalar *x,*d;
1558fe8eb27SJed Brown     PetscScalar       *y;
1568fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
1578fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1588fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1598fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1608fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
1618fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1628fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1638fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
164880538faSHong Zhang   } else if (PetscAbsScalar(shell->vshift) != 0) {
1658fe8eb27SJed Brown     ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr);
166027c5fb5SJed Brown   } else if (shell->vscale != 1.0) {
167027c5fb5SJed Brown     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
1688fe8eb27SJed Brown   }
1698fe8eb27SJed Brown   PetscFunctionReturn(0);
1708fe8eb27SJed Brown }
171e51e0e81SBarry Smith 
172711e205bSSatish Balay #undef __FUNCT__
173711e205bSSatish Balay #define __FUNCT__ "MatShellGetContext"
1749d225801SJed Brown /*@
175a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
176b4fd4287SBarry Smith 
177c7fcc2eaSBarry Smith     Not Collective
178c7fcc2eaSBarry Smith 
179b4fd4287SBarry Smith     Input Parameter:
180b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
181b4fd4287SBarry Smith 
182b4fd4287SBarry Smith     Output Parameter:
183b4fd4287SBarry Smith .   ctx - the user provided context
184b4fd4287SBarry Smith 
18515091d37SBarry Smith     Level: advanced
18615091d37SBarry Smith 
187a62d957aSLois Curfman McInnes     Notes:
188a62d957aSLois Curfman McInnes     This routine is intended for use within various shell matrix routines,
189a62d957aSLois Curfman McInnes     as set with MatShellSetOperation().
190a62d957aSLois Curfman McInnes 
191a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context
192a62d957aSLois Curfman McInnes 
193ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
1940bc0a719SBarry Smith @*/
1958fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
196b4fd4287SBarry Smith {
197dfbe8321SBarry Smith   PetscErrorCode ierr;
198ace3abfcSBarry Smith   PetscBool      flg;
199273d9f13SBarry Smith 
2003a40ed3dSBarry Smith   PetscFunctionBegin;
2010700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2024482741eSBarry Smith   PetscValidPointer(ctx,2);
203251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
204940183f3SBarry Smith   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
205940183f3SBarry Smith   else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
2063a40ed3dSBarry Smith   PetscFunctionReturn(0);
207b4fd4287SBarry Smith }
208b4fd4287SBarry Smith 
209711e205bSSatish Balay #undef __FUNCT__
210711e205bSSatish Balay #define __FUNCT__ "MatDestroy_Shell"
211dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
212e51e0e81SBarry Smith {
213dfbe8321SBarry Smith   PetscErrorCode ierr;
214bf0cc555SLisandro Dalcin   Mat_Shell      *shell = (Mat_Shell*)mat->data;
215ed3cc1f0SBarry Smith 
2163a40ed3dSBarry Smith   PetscFunctionBegin;
217bf0cc555SLisandro Dalcin   if (shell->destroy) {
218bf0cc555SLisandro Dalcin     ierr = (*shell->destroy)(mat);CHKERRQ(ierr);
219bf0cc555SLisandro Dalcin   }
2208fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr);
2218fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr);
2228fe8eb27SJed Brown   ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr);
2238fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
2248fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
225*5edf6516SJed Brown   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
226*5edf6516SJed Brown   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
227bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
2283a40ed3dSBarry Smith   PetscFunctionReturn(0);
229e51e0e81SBarry Smith }
230e51e0e81SBarry Smith 
231711e205bSSatish Balay #undef __FUNCT__
232ef66eb69SBarry Smith #define __FUNCT__ "MatMult_Shell"
233dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
234ef66eb69SBarry Smith {
235ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)A->data;
236dfbe8321SBarry Smith   PetscErrorCode ierr;
23725578ef6SJed Brown   Vec            xx;
238ef66eb69SBarry Smith 
239ef66eb69SBarry Smith   PetscFunctionBegin;
2408fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
2418fe8eb27SJed Brown   ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr);
2428fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2438fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
244ef66eb69SBarry Smith   PetscFunctionReturn(0);
245ef66eb69SBarry Smith }
246ef66eb69SBarry Smith 
247ef66eb69SBarry Smith #undef __FUNCT__
248*5edf6516SJed Brown #define __FUNCT__ "MatMultAdd_Shell"
249*5edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
250*5edf6516SJed Brown {
251*5edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
252*5edf6516SJed Brown   PetscErrorCode ierr;
253*5edf6516SJed Brown 
254*5edf6516SJed Brown   PetscFunctionBegin;
255*5edf6516SJed Brown   if (y == z) {
256*5edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
257*5edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
258*5edf6516SJed Brown     ierr = VecWAXPY(z,1.0,shell->right_add_work,y);CHKERRQ(ierr);
259*5edf6516SJed Brown   } else {
260*5edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
261*5edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
262*5edf6516SJed Brown   }
263*5edf6516SJed Brown   PetscFunctionReturn(0);
264*5edf6516SJed Brown }
265*5edf6516SJed Brown 
266*5edf6516SJed Brown #undef __FUNCT__
26718be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell"
26818be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
26918be62a5SSatish Balay {
27018be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
27118be62a5SSatish Balay   PetscErrorCode ierr;
27225578ef6SJed Brown   Vec            xx;
27318be62a5SSatish Balay 
27418be62a5SSatish Balay   PetscFunctionBegin;
2758fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
2768fe8eb27SJed Brown   ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr);
2778fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2788fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
27918be62a5SSatish Balay   PetscFunctionReturn(0);
28018be62a5SSatish Balay }
28118be62a5SSatish Balay 
28218be62a5SSatish Balay #undef __FUNCT__
283*5edf6516SJed Brown #define __FUNCT__ "MatMultTransposeAdd_Shell"
284*5edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
285*5edf6516SJed Brown {
286*5edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
287*5edf6516SJed Brown   PetscErrorCode ierr;
288*5edf6516SJed Brown 
289*5edf6516SJed Brown   PetscFunctionBegin;
290*5edf6516SJed Brown   if (y == z) {
291*5edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
292*5edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
293*5edf6516SJed Brown     ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr);
294*5edf6516SJed Brown   } else {
295*5edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
296*5edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
297*5edf6516SJed Brown   }
298*5edf6516SJed Brown   PetscFunctionReturn(0);
299*5edf6516SJed Brown }
300*5edf6516SJed Brown 
301*5edf6516SJed Brown #undef __FUNCT__
30218be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell"
30318be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
30418be62a5SSatish Balay {
30518be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
30618be62a5SSatish Balay   PetscErrorCode ierr;
30718be62a5SSatish Balay 
30818be62a5SSatish Balay   PetscFunctionBegin;
30918be62a5SSatish Balay   ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr);
31018be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
3118fe8eb27SJed Brown   if (shell->dshift) {
3128fe8eb27SJed Brown     ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr);
3138fe8eb27SJed Brown   } else {
31418be62a5SSatish Balay     ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
31518be62a5SSatish Balay   }
3168fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
3178fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
31818be62a5SSatish Balay   PetscFunctionReturn(0);
31918be62a5SSatish Balay }
32018be62a5SSatish Balay 
32118be62a5SSatish Balay #undef __FUNCT__
322ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell"
323f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
324ef66eb69SBarry Smith {
325ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
3268fe8eb27SJed Brown   PetscErrorCode ierr;
327b24ad042SBarry Smith 
328ef66eb69SBarry Smith   PetscFunctionBegin;
3298fe8eb27SJed Brown   if (shell->left || shell->right || shell->dshift) {
3308fe8eb27SJed Brown     if (!shell->dshift) {
3318fe8eb27SJed Brown       if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);}
3328fe8eb27SJed Brown       shell->dshift = shell->dshift_owned;
3338fe8eb27SJed Brown       ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr);
3348fe8eb27SJed Brown     } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
3358fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
3368fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
3378fe8eb27SJed Brown   } else shell->vshift += a;
3388fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
339ef66eb69SBarry Smith   PetscFunctionReturn(0);
340ef66eb69SBarry Smith }
341ef66eb69SBarry Smith 
342ef66eb69SBarry Smith #undef __FUNCT__
343ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell"
344f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
345ef66eb69SBarry Smith {
346ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3478fe8eb27SJed Brown   PetscErrorCode ierr;
348b24ad042SBarry Smith 
349ef66eb69SBarry Smith   PetscFunctionBegin;
350f4df32b1SMatthew Knepley   shell->vscale *= a;
3518fe8eb27SJed Brown   if (shell->dshift) {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
3528fe8eb27SJed Brown   else shell->vshift *= a;
3538fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
3548fe8eb27SJed Brown   PetscFunctionReturn(0);
35518be62a5SSatish Balay }
3568fe8eb27SJed Brown 
3578fe8eb27SJed Brown #undef __FUNCT__
3588fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell"
3598fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
3608fe8eb27SJed Brown {
3618fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)Y->data;
3628fe8eb27SJed Brown   PetscErrorCode ierr;
3638fe8eb27SJed Brown 
3648fe8eb27SJed Brown   PetscFunctionBegin;
3658fe8eb27SJed Brown   if (left) {
3668fe8eb27SJed Brown     if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);}
3678fe8eb27SJed Brown     if (shell->left) {ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);}
3688fe8eb27SJed Brown     else {
3698fe8eb27SJed Brown       shell->left = shell->left_owned;
3708fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
37118be62a5SSatish Balay     }
372ef66eb69SBarry Smith   }
3738fe8eb27SJed Brown   if (right) {
3748fe8eb27SJed Brown     if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);}
3758fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);}
3768fe8eb27SJed Brown     else {
3778fe8eb27SJed Brown       shell->right = shell->right_owned;
3788fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
3798fe8eb27SJed Brown     }
3808fe8eb27SJed Brown   }
3818fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
382ef66eb69SBarry Smith   PetscFunctionReturn(0);
383ef66eb69SBarry Smith }
384ef66eb69SBarry Smith 
385ef66eb69SBarry Smith #undef __FUNCT__
386ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell"
387dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
388ef66eb69SBarry Smith {
389ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
390ef66eb69SBarry Smith 
391ef66eb69SBarry Smith   PetscFunctionBegin;
3928fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
393ef66eb69SBarry Smith     shell->vshift      = 0.0;
394ef66eb69SBarry Smith     shell->vscale      = 1.0;
3958fe8eb27SJed Brown     shell->dshift      = PETSC_NULL;
3968fe8eb27SJed Brown     shell->left        = PETSC_NULL;
3978fe8eb27SJed Brown     shell->right       = PETSC_NULL;
3987fb7d96aSJed Brown     if (shell->mult) {
399ef66eb69SBarry Smith       Y->ops->mult = shell->mult;
4007fb7d96aSJed Brown       shell->mult  = PETSC_NULL;
4017fb7d96aSJed Brown     }
4027fb7d96aSJed Brown     if (shell->multtranspose) {
40318be62a5SSatish Balay       Y->ops->multtranspose = shell->multtranspose;
4047fb7d96aSJed Brown       shell->multtranspose  = PETSC_NULL;
4057fb7d96aSJed Brown     }
4067fb7d96aSJed Brown     if (shell->getdiagonal) {
40718be62a5SSatish Balay       Y->ops->getdiagonal = shell->getdiagonal;
4087fb7d96aSJed Brown       shell->getdiagonal  = PETSC_NULL;
4097fb7d96aSJed Brown     }
4107fb7d96aSJed Brown     shell->usingscaled = PETSC_FALSE;
411ef66eb69SBarry Smith   }
412ef66eb69SBarry Smith   PetscFunctionReturn(0);
413ef66eb69SBarry Smith }
414ef66eb69SBarry Smith 
41519fd82e9SBarry Smith extern PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
416b951964fSBarry Smith 
41709dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
41820563c6bSBarry Smith        0,
41920563c6bSBarry Smith        0,
42020563c6bSBarry Smith        0,
42197304618SKris Buschelman /* 4*/ 0,
42220563c6bSBarry Smith        0,
423b951964fSBarry Smith        0,
424b951964fSBarry Smith        0,
425b951964fSBarry Smith        0,
426b951964fSBarry Smith        0,
42797304618SKris Buschelman /*10*/ 0,
428b951964fSBarry Smith        0,
429b951964fSBarry Smith        0,
430b951964fSBarry Smith        0,
431b951964fSBarry Smith        0,
43297304618SKris Buschelman /*15*/ 0,
433b951964fSBarry Smith        0,
434b951964fSBarry Smith        0,
4358fe8eb27SJed Brown        MatDiagonalScale_Shell,
436b951964fSBarry Smith        0,
43797304618SKris Buschelman /*20*/ 0,
438ef66eb69SBarry Smith        MatAssemblyEnd_Shell,
439b951964fSBarry Smith        0,
440b951964fSBarry Smith        0,
441d519adbfSMatthew Knepley /*24*/ 0,
442b951964fSBarry Smith        0,
443b951964fSBarry Smith        0,
444b951964fSBarry Smith        0,
445b951964fSBarry Smith        0,
446d519adbfSMatthew Knepley /*29*/ 0,
447b951964fSBarry Smith        0,
448273d9f13SBarry Smith        0,
449b951964fSBarry Smith        0,
450b951964fSBarry Smith        0,
451d519adbfSMatthew Knepley /*34*/ 0,
452b951964fSBarry Smith        0,
453b951964fSBarry Smith        0,
45409dc0095SBarry Smith        0,
45509dc0095SBarry Smith        0,
456d519adbfSMatthew Knepley /*39*/ 0,
45709dc0095SBarry Smith        0,
45809dc0095SBarry Smith        0,
45909dc0095SBarry Smith        0,
46009dc0095SBarry Smith        0,
461d519adbfSMatthew Knepley /*44*/ 0,
462ef66eb69SBarry Smith        MatScale_Shell,
463ef66eb69SBarry Smith        MatShift_Shell,
46409dc0095SBarry Smith        0,
46509dc0095SBarry Smith        0,
466f73d5cc4SBarry Smith /*49*/ 0,
46709dc0095SBarry Smith        0,
46809dc0095SBarry Smith        0,
46909dc0095SBarry Smith        0,
47009dc0095SBarry Smith        0,
471d519adbfSMatthew Knepley /*54*/ 0,
47209dc0095SBarry Smith        0,
47309dc0095SBarry Smith        0,
47409dc0095SBarry Smith        0,
47509dc0095SBarry Smith        0,
476d519adbfSMatthew Knepley /*59*/ 0,
477b9b97703SBarry Smith        MatDestroy_Shell,
47809dc0095SBarry Smith        0,
479357abbc8SBarry Smith        0,
480273d9f13SBarry Smith        0,
481d519adbfSMatthew Knepley /*64*/ 0,
482273d9f13SBarry Smith        0,
483273d9f13SBarry Smith        0,
484273d9f13SBarry Smith        0,
485273d9f13SBarry Smith        0,
486d519adbfSMatthew Knepley /*69*/ 0,
487273d9f13SBarry Smith        0,
488c87e5d42SMatthew Knepley        MatConvert_Shell,
489273d9f13SBarry Smith        0,
49097304618SKris Buschelman        0,
491d519adbfSMatthew Knepley /*74*/ 0,
49297304618SKris Buschelman        0,
49397304618SKris Buschelman        0,
49497304618SKris Buschelman        0,
49597304618SKris Buschelman        0,
496d519adbfSMatthew Knepley /*79*/ 0,
49797304618SKris Buschelman        0,
49897304618SKris Buschelman        0,
49997304618SKris Buschelman        0,
500865e5f61SKris Buschelman        0,
501d519adbfSMatthew Knepley /*84*/ 0,
502865e5f61SKris Buschelman        0,
503865e5f61SKris Buschelman        0,
504865e5f61SKris Buschelman        0,
505865e5f61SKris Buschelman        0,
506d519adbfSMatthew Knepley /*89*/ 0,
507865e5f61SKris Buschelman        0,
508865e5f61SKris Buschelman        0,
509865e5f61SKris Buschelman        0,
510865e5f61SKris Buschelman        0,
511d519adbfSMatthew Knepley /*94*/ 0,
512865e5f61SKris Buschelman        0,
513865e5f61SKris Buschelman        0,
514865e5f61SKris Buschelman        0};
515273d9f13SBarry Smith 
5160bad9183SKris Buschelman /*MC
517fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
5180bad9183SKris Buschelman 
5190bad9183SKris Buschelman   Level: advanced
5200bad9183SKris Buschelman 
5210bad9183SKris Buschelman .seealso: MatCreateShell
5220bad9183SKris Buschelman M*/
5230bad9183SKris Buschelman 
524273d9f13SBarry Smith EXTERN_C_BEGIN
525711e205bSSatish Balay #undef __FUNCT__
526711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell"
5277087cfbeSBarry Smith PetscErrorCode  MatCreate_Shell(Mat A)
528273d9f13SBarry Smith {
529273d9f13SBarry Smith   Mat_Shell      *b;
530dfbe8321SBarry Smith   PetscErrorCode ierr;
531273d9f13SBarry Smith 
532273d9f13SBarry Smith   PetscFunctionBegin;
533273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
534273d9f13SBarry Smith 
53538f2d2fdSLisandro Dalcin   ierr = PetscNewLog(A,Mat_Shell,&b);CHKERRQ(ierr);
536273d9f13SBarry Smith   A->data = (void*)b;
537273d9f13SBarry Smith 
53826283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
53926283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
540273d9f13SBarry Smith 
541273d9f13SBarry Smith   b->ctx           = 0;
542ef66eb69SBarry Smith   b->vshift        = 0.0;
543ef66eb69SBarry Smith   b->vscale        = 1.0;
544ef66eb69SBarry Smith   b->mult          = 0;
54518be62a5SSatish Balay   b->multtranspose = 0;
54618be62a5SSatish Balay   b->getdiagonal   = 0;
547273d9f13SBarry Smith   A->assembled     = PETSC_TRUE;
548210f0121SBarry Smith   A->preallocated  = PETSC_FALSE;
54917667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
550273d9f13SBarry Smith   PetscFunctionReturn(0);
551273d9f13SBarry Smith }
552273d9f13SBarry Smith EXTERN_C_END
553e51e0e81SBarry Smith 
554711e205bSSatish Balay #undef __FUNCT__
555711e205bSSatish Balay #define __FUNCT__ "MatCreateShell"
5564b828684SBarry Smith /*@C
557052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
558ff756334SLois Curfman McInnes    private data storage format.
559e51e0e81SBarry Smith 
560c7fcc2eaSBarry Smith   Collective on MPI_Comm
561c7fcc2eaSBarry Smith 
562e51e0e81SBarry Smith    Input Parameters:
563c7fcc2eaSBarry Smith +  comm - MPI communicator
564c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
565c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
566c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
567c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
568c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
569e51e0e81SBarry Smith 
570ff756334SLois Curfman McInnes    Output Parameter:
57144cd7ae7SLois Curfman McInnes .  A - the matrix
572e51e0e81SBarry Smith 
573ff2fd236SBarry Smith    Level: advanced
574ff2fd236SBarry Smith 
575f39d1f56SLois Curfman McInnes   Usage:
5767b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
577f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
578c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
579f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
580f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
581f39d1f56SLois Curfman McInnes 
582ff756334SLois Curfman McInnes    Notes:
583ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
584ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
585ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
586e51e0e81SBarry Smith 
587f38a8d56SBarry Smith    Fortran Notes: The context can only be an integer or a PetscObject
588f38a8d56SBarry Smith       unfortunately it cannot be a Fortran array or derived type.
589f38a8d56SBarry Smith 
590f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
591f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
592645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
593645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
594645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
595645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
596645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
597645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
598645985a0SLois Curfman McInnes    For example,
599f39d1f56SLois Curfman McInnes 
600f39d1f56SLois Curfman McInnes $
601f39d1f56SLois Curfman McInnes $     Vec x, y
6027b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
603645985a0SLois Curfman McInnes $     Mat A
604f39d1f56SLois Curfman McInnes $
605c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
606c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
607f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
608c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
609c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
610c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
611645985a0SLois Curfman McInnes $     MatMult(A,x,y);
612645985a0SLois Curfman McInnes $     MatDestroy(A);
613f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
614645985a0SLois Curfman McInnes $
615e51e0e81SBarry Smith 
6160b627109SLois Curfman McInnes .keywords: matrix, shell, create
6170b627109SLois Curfman McInnes 
618ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
619e51e0e81SBarry Smith @*/
6207087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
621e51e0e81SBarry Smith {
622dfbe8321SBarry Smith   PetscErrorCode ierr;
623ed3cc1f0SBarry Smith 
6243a40ed3dSBarry Smith   PetscFunctionBegin;
625f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
626f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
627273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
628273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
6290fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
630273d9f13SBarry Smith   PetscFunctionReturn(0);
631c7fcc2eaSBarry Smith }
632c7fcc2eaSBarry Smith 
633711e205bSSatish Balay #undef __FUNCT__
634711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext"
635c6866cfdSSatish Balay /*@
636273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
637c7fcc2eaSBarry Smith 
6383f9fe445SBarry Smith    Logically Collective on Mat
639c7fcc2eaSBarry Smith 
640273d9f13SBarry Smith     Input Parameters:
641273d9f13SBarry Smith +   mat - the shell matrix
642273d9f13SBarry Smith -   ctx - the context
643273d9f13SBarry Smith 
644273d9f13SBarry Smith    Level: advanced
645273d9f13SBarry Smith 
646f38a8d56SBarry Smith    Fortran Notes: The context can only be an integer or a PetscObject
647f38a8d56SBarry Smith       unfortunately it cannot be a Fortran array or derived type.
648273d9f13SBarry Smith 
649273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
6500bc0a719SBarry Smith @*/
6517087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
652273d9f13SBarry Smith {
653273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
654dfbe8321SBarry Smith   PetscErrorCode ierr;
655ace3abfcSBarry Smith   PetscBool      flg;
656273d9f13SBarry Smith 
657273d9f13SBarry Smith   PetscFunctionBegin;
6580700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
659251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
660273d9f13SBarry Smith   if (flg) {
661273d9f13SBarry Smith     shell->ctx = ctx;
662940183f3SBarry Smith   } else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
6633a40ed3dSBarry Smith   PetscFunctionReturn(0);
664e51e0e81SBarry Smith }
665e51e0e81SBarry Smith 
666711e205bSSatish Balay #undef __FUNCT__
667711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation"
668c16cb8f2SBarry Smith /*@C
6693a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
6703a3eedf2SBarry Smith                            a shell matrix.
671e51e0e81SBarry Smith 
6723f9fe445SBarry Smith    Logically Collective on Mat
673fee21e36SBarry Smith 
674c7fcc2eaSBarry Smith     Input Parameters:
675c7fcc2eaSBarry Smith +   mat - the shell matrix
676c7fcc2eaSBarry Smith .   op - the name of the operation
677c7fcc2eaSBarry Smith -   f - the function that provides the operation.
678c7fcc2eaSBarry Smith 
67915091d37SBarry Smith    Level: advanced
68015091d37SBarry Smith 
681fae171e0SBarry Smith     Usage:
682e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
683f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
684c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
6850b627109SLois Curfman McInnes 
686a62d957aSLois Curfman McInnes     Notes:
687e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
6881c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
689a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
6901c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
691a62d957aSLois Curfman McInnes 
69225296bd5SBarry Smith     All user-provided functions (execept for MATOP_DESTROY) should have the same calling
693deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
694deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
695deebb3c3SLois Curfman McInnes     routines, e.g.,
696a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
697a62d957aSLois Curfman McInnes 
6984aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
6994aa34b0aSBarry Smith     nonzero on failure.
7004aa34b0aSBarry Smith 
701a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
702a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
703a62d957aSLois Curfman McInnes     set by MatCreateShell().
704a62d957aSLois Curfman McInnes 
705501d9185SBarry Smith     Fortran Notes: For MatGetVecs() the user code should check if the input left or right matrix is -1 and in that case not
706501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
707501d9185SBarry Smith 
708a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
709a62d957aSLois Curfman McInnes 
710ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
711e51e0e81SBarry Smith @*/
7127087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
713e51e0e81SBarry Smith {
714dfbe8321SBarry Smith   PetscErrorCode ierr;
715ace3abfcSBarry Smith   PetscBool      flg;
716273d9f13SBarry Smith 
7173a40ed3dSBarry Smith   PetscFunctionBegin;
7180700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
719*5edf6516SJed Brown   switch (op) {
720*5edf6516SJed Brown   case MATOP_DESTROY:
721251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
722273d9f13SBarry Smith     if (flg) {
723a62d957aSLois Curfman McInnes        Mat_Shell *shell = (Mat_Shell*)mat->data;
7246849ba73SBarry Smith        shell->destroy                 = (PetscErrorCode (*)(Mat)) f;
7256849ba73SBarry Smith     } else mat->ops->destroy          = (PetscErrorCode (*)(Mat)) f;
726*5edf6516SJed Brown     break;
727*5edf6516SJed Brown   case MATOP_VIEW:
728*5edf6516SJed Brown     mat->ops->view  = (PetscErrorCode (*)(Mat,PetscViewer)) f;
729*5edf6516SJed Brown     break;
730*5edf6516SJed Brown   case MATOP_MULT:
731*5edf6516SJed Brown     mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec)) f;
732*5edf6516SJed Brown     if (!mat->ops->multadd) mat->ops->multadd = MatMultAdd_Shell;
733*5edf6516SJed Brown     break;
734*5edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
735*5edf6516SJed Brown     mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec)) f;
736*5edf6516SJed Brown     if (!mat->ops->multtransposeadd) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell;
737*5edf6516SJed Brown     break;
738*5edf6516SJed Brown   default:
739*5edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
740a62d957aSLois Curfman McInnes   }
7413a40ed3dSBarry Smith   PetscFunctionReturn(0);
742e51e0e81SBarry Smith }
743f0479e8cSBarry Smith 
744711e205bSSatish Balay #undef __FUNCT__
745711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation"
746d4bb536fSBarry Smith /*@C
747d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
748d4bb536fSBarry Smith 
749c7fcc2eaSBarry Smith     Not Collective
750c7fcc2eaSBarry Smith 
751d4bb536fSBarry Smith     Input Parameters:
752c7fcc2eaSBarry Smith +   mat - the shell matrix
753c7fcc2eaSBarry Smith -   op - the name of the operation
754d4bb536fSBarry Smith 
755d4bb536fSBarry Smith     Output Parameter:
756d4bb536fSBarry Smith .   f - the function that provides the operation.
757d4bb536fSBarry Smith 
75815091d37SBarry Smith     Level: advanced
75915091d37SBarry Smith 
760d4bb536fSBarry Smith     Notes:
761e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
762d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
763d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
764d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
765d4bb536fSBarry Smith 
766d4bb536fSBarry Smith     All user-provided functions have the same calling
767d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
768d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
769d4bb536fSBarry Smith     routines, e.g.,
770d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
771d4bb536fSBarry Smith 
772d4bb536fSBarry Smith     Within each user-defined routine, the user should call
773d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
774d4bb536fSBarry Smith     set by MatCreateShell().
775d4bb536fSBarry Smith 
776d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
777d4bb536fSBarry Smith 
778ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
779d4bb536fSBarry Smith @*/
7807087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
781d4bb536fSBarry Smith {
782dfbe8321SBarry Smith   PetscErrorCode ierr;
783ace3abfcSBarry Smith   PetscBool      flg;
784273d9f13SBarry Smith 
7853a40ed3dSBarry Smith   PetscFunctionBegin;
7860700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
787d4bb536fSBarry Smith   if (op == MATOP_DESTROY) {
788251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
789273d9f13SBarry Smith     if (flg) {
790d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
791c134de8dSSatish Balay       *f = (void(*)(void))shell->destroy;
792c7fcc2eaSBarry Smith     } else {
793c134de8dSSatish Balay       *f = (void(*)(void))mat->ops->destroy;
794d4bb536fSBarry Smith     }
795c7fcc2eaSBarry Smith   } else if (op == MATOP_VIEW) {
796c134de8dSSatish Balay     *f = (void(*)(void))mat->ops->view;
797c7fcc2eaSBarry Smith   } else {
798c134de8dSSatish Balay     *f = (((void(**)(void))mat->ops)[op]);
799d4bb536fSBarry Smith   }
800d4bb536fSBarry Smith 
8013a40ed3dSBarry Smith   PetscFunctionReturn(0);
802d4bb536fSBarry Smith }
803d4bb536fSBarry Smith 
804