xref: /petsc/src/mat/impls/shell/shell.c (revision 251f4c6745e69799cb0300fd80598c58f96eb452)
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;
218fe8eb27SJed Brown   PetscBool      usingscaled;
2220563c6bSBarry Smith   void           *ctx;
2388cf3e7dSBarry Smith } Mat_Shell;
248fe8eb27SJed Brown /*
258fe8eb27SJed Brown  The most general expression for the matrix is
268fe8eb27SJed Brown 
278fe8eb27SJed Brown  S = L (a A + B) R
288fe8eb27SJed Brown 
298fe8eb27SJed Brown  where
308fe8eb27SJed Brown  A is the matrix defined by the user's function
318fe8eb27SJed Brown  a is a scalar multiple
328fe8eb27SJed Brown  L is left scaling
338fe8eb27SJed Brown  R is right scaling
348fe8eb27SJed Brown  B is a diagonal shift defined by
358fe8eb27SJed Brown    diag(dshift) if the vector dshift is non-NULL
368fe8eb27SJed Brown    vscale*identity otherwise
378fe8eb27SJed Brown 
388fe8eb27SJed Brown  The following identities apply:
398fe8eb27SJed Brown 
408fe8eb27SJed Brown  Scale by c:
418fe8eb27SJed Brown   c [L (a A + B) R] = L [(a c) A + c B] R
428fe8eb27SJed Brown 
438fe8eb27SJed Brown  Shift by c:
448fe8eb27SJed Brown   [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R
458fe8eb27SJed Brown 
468fe8eb27SJed Brown  Diagonal scaling is achieved by simply multiplying with existing L and R vectors
478fe8eb27SJed Brown 
488fe8eb27SJed Brown  In the data structure:
498fe8eb27SJed Brown 
508fe8eb27SJed Brown  vscale=1.0  means no special scaling will be applied
518fe8eb27SJed Brown  dshift=NULL means a constant diagonal shift (fall back to vshift)
528fe8eb27SJed Brown  vshift=0.0  means no constant diagonal shift, note that vshift is only used if dshift is NULL
538fe8eb27SJed Brown */
548fe8eb27SJed Brown 
558fe8eb27SJed Brown static PetscErrorCode MatMult_Shell(Mat,Vec,Vec);
568fe8eb27SJed Brown static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec);
578fe8eb27SJed Brown static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec);
588fe8eb27SJed Brown 
598fe8eb27SJed Brown #undef __FUNCT__
608fe8eb27SJed Brown #define __FUNCT__ "MatShellUseScaledMethods"
618fe8eb27SJed Brown static PetscErrorCode MatShellUseScaledMethods(Mat Y)
628fe8eb27SJed Brown {
638fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)Y->data;
648fe8eb27SJed Brown 
658fe8eb27SJed Brown   PetscFunctionBegin;
668fe8eb27SJed Brown   if (shell->usingscaled) PetscFunctionReturn(0);
678fe8eb27SJed Brown   shell->mult  = Y->ops->mult;
688fe8eb27SJed Brown   Y->ops->mult = MatMult_Shell;
698fe8eb27SJed Brown   if (Y->ops->multtranspose) {
708fe8eb27SJed Brown     shell->multtranspose  = Y->ops->multtranspose;
718fe8eb27SJed Brown     Y->ops->multtranspose = MatMultTranspose_Shell;
728fe8eb27SJed Brown   }
738fe8eb27SJed Brown   if (Y->ops->getdiagonal) {
748fe8eb27SJed Brown     shell->getdiagonal  = Y->ops->getdiagonal;
758fe8eb27SJed Brown     Y->ops->getdiagonal = MatGetDiagonal_Shell;
768fe8eb27SJed Brown   }
778fe8eb27SJed Brown   shell->usingscaled = PETSC_TRUE;
788fe8eb27SJed Brown   PetscFunctionReturn(0);
798fe8eb27SJed Brown }
808fe8eb27SJed Brown 
818fe8eb27SJed Brown #undef __FUNCT__
828fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleLeft"
838fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
848fe8eb27SJed Brown {
858fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
868fe8eb27SJed Brown   PetscErrorCode ierr;
878fe8eb27SJed Brown 
888fe8eb27SJed Brown   PetscFunctionBegin;
89d63f877fSJed Brown   *xx = PETSC_NULL;
908fe8eb27SJed Brown   if (!shell->left) {
918fe8eb27SJed Brown     *xx = x;
928fe8eb27SJed Brown   } else {
938fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
948fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
958fe8eb27SJed Brown     *xx = shell->left_work;
968fe8eb27SJed Brown   }
978fe8eb27SJed Brown   PetscFunctionReturn(0);
988fe8eb27SJed Brown }
998fe8eb27SJed Brown 
1008fe8eb27SJed Brown #undef __FUNCT__
1018fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleRight"
1028fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
1038fe8eb27SJed Brown {
1048fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1058fe8eb27SJed Brown   PetscErrorCode ierr;
1068fe8eb27SJed Brown 
1078fe8eb27SJed Brown   PetscFunctionBegin;
108d63f877fSJed Brown   *xx = PETSC_NULL;
1098fe8eb27SJed Brown   if (!shell->right) {
1108fe8eb27SJed Brown     *xx = x;
1118fe8eb27SJed Brown   } else {
1128fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
1138fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
1148fe8eb27SJed Brown     *xx = shell->right_work;
1158fe8eb27SJed Brown   }
1168fe8eb27SJed Brown   PetscFunctionReturn(0);
1178fe8eb27SJed Brown }
1188fe8eb27SJed Brown 
1198fe8eb27SJed Brown #undef __FUNCT__
1208fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleLeft"
1218fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
1228fe8eb27SJed Brown {
1238fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1248fe8eb27SJed Brown   PetscErrorCode ierr;
1258fe8eb27SJed Brown 
1268fe8eb27SJed Brown   PetscFunctionBegin;
1278fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
1288fe8eb27SJed Brown   PetscFunctionReturn(0);
1298fe8eb27SJed Brown }
1308fe8eb27SJed Brown 
1318fe8eb27SJed Brown #undef __FUNCT__
1328fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleRight"
1338fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
1348fe8eb27SJed Brown {
1358fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1368fe8eb27SJed Brown   PetscErrorCode ierr;
1378fe8eb27SJed Brown 
1388fe8eb27SJed Brown   PetscFunctionBegin;
1398fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
1408fe8eb27SJed Brown   PetscFunctionReturn(0);
1418fe8eb27SJed Brown }
1428fe8eb27SJed Brown 
1438fe8eb27SJed Brown #undef __FUNCT__
1448fe8eb27SJed Brown #define __FUNCT__ "MatShellShiftAndScale"
1458fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
1468fe8eb27SJed Brown {
1478fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1488fe8eb27SJed Brown   PetscErrorCode ierr;
1498fe8eb27SJed Brown 
1508fe8eb27SJed Brown   PetscFunctionBegin;
1518fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
1528fe8eb27SJed Brown     PetscInt          i,m;
1538fe8eb27SJed Brown     const PetscScalar *x,*d;
1548fe8eb27SJed Brown     PetscScalar       *y;
1558fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
1568fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1578fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1588fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1598fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
1608fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1618fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1628fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
163880538faSHong Zhang   } else if (PetscAbsScalar(shell->vshift) != 0) {
1648fe8eb27SJed Brown     ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr);
165027c5fb5SJed Brown   } else if (shell->vscale != 1.0) {
166027c5fb5SJed Brown     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
1678fe8eb27SJed Brown   }
1688fe8eb27SJed Brown   PetscFunctionReturn(0);
1698fe8eb27SJed Brown }
170e51e0e81SBarry Smith 
171711e205bSSatish Balay #undef __FUNCT__
172711e205bSSatish Balay #define __FUNCT__ "MatShellGetContext"
1739d225801SJed Brown /*@
174a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
175b4fd4287SBarry Smith 
176c7fcc2eaSBarry Smith     Not Collective
177c7fcc2eaSBarry Smith 
178b4fd4287SBarry Smith     Input Parameter:
179b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
180b4fd4287SBarry Smith 
181b4fd4287SBarry Smith     Output Parameter:
182b4fd4287SBarry Smith .   ctx - the user provided context
183b4fd4287SBarry Smith 
18415091d37SBarry Smith     Level: advanced
18515091d37SBarry Smith 
186a62d957aSLois Curfman McInnes     Notes:
187a62d957aSLois Curfman McInnes     This routine is intended for use within various shell matrix routines,
188a62d957aSLois Curfman McInnes     as set with MatShellSetOperation().
189a62d957aSLois Curfman McInnes 
190a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context
191a62d957aSLois Curfman McInnes 
192ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
1930bc0a719SBarry Smith @*/
1948fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
195b4fd4287SBarry Smith {
196dfbe8321SBarry Smith   PetscErrorCode ierr;
197ace3abfcSBarry Smith   PetscBool      flg;
198273d9f13SBarry Smith 
1993a40ed3dSBarry Smith   PetscFunctionBegin;
2000700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2014482741eSBarry Smith   PetscValidPointer(ctx,2);
202*251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
203940183f3SBarry Smith   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
204940183f3SBarry Smith   else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
2053a40ed3dSBarry Smith   PetscFunctionReturn(0);
206b4fd4287SBarry Smith }
207b4fd4287SBarry Smith 
208711e205bSSatish Balay #undef __FUNCT__
209711e205bSSatish Balay #define __FUNCT__ "MatDestroy_Shell"
210dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
211e51e0e81SBarry Smith {
212dfbe8321SBarry Smith   PetscErrorCode ierr;
213bf0cc555SLisandro Dalcin   Mat_Shell      *shell = (Mat_Shell*)mat->data;
214ed3cc1f0SBarry Smith 
2153a40ed3dSBarry Smith   PetscFunctionBegin;
216bf0cc555SLisandro Dalcin   if (shell->destroy) {
217bf0cc555SLisandro Dalcin     ierr = (*shell->destroy)(mat);CHKERRQ(ierr);
218bf0cc555SLisandro Dalcin   }
2198fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr);
2208fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr);
2218fe8eb27SJed Brown   ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr);
2228fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
2238fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
224bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
2253a40ed3dSBarry Smith   PetscFunctionReturn(0);
226e51e0e81SBarry Smith }
227e51e0e81SBarry Smith 
228711e205bSSatish Balay #undef __FUNCT__
229ef66eb69SBarry Smith #define __FUNCT__ "MatMult_Shell"
230dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
231ef66eb69SBarry Smith {
232ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)A->data;
233dfbe8321SBarry Smith   PetscErrorCode ierr;
23425578ef6SJed Brown   Vec            xx;
235ef66eb69SBarry Smith 
236ef66eb69SBarry Smith   PetscFunctionBegin;
2378fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
2388fe8eb27SJed Brown   ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr);
2398fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2408fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
241ef66eb69SBarry Smith   PetscFunctionReturn(0);
242ef66eb69SBarry Smith }
243ef66eb69SBarry Smith 
244ef66eb69SBarry Smith #undef __FUNCT__
24518be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell"
24618be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
24718be62a5SSatish Balay {
24818be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
24918be62a5SSatish Balay   PetscErrorCode ierr;
25025578ef6SJed Brown   Vec            xx;
25118be62a5SSatish Balay 
25218be62a5SSatish Balay   PetscFunctionBegin;
2538fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
2548fe8eb27SJed Brown   ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr);
2558fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2568fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
25718be62a5SSatish Balay   PetscFunctionReturn(0);
25818be62a5SSatish Balay }
25918be62a5SSatish Balay 
26018be62a5SSatish Balay #undef __FUNCT__
26118be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell"
26218be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
26318be62a5SSatish Balay {
26418be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
26518be62a5SSatish Balay   PetscErrorCode ierr;
26618be62a5SSatish Balay 
26718be62a5SSatish Balay   PetscFunctionBegin;
26818be62a5SSatish Balay   ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr);
26918be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
2708fe8eb27SJed Brown   if (shell->dshift) {
2718fe8eb27SJed Brown     ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr);
2728fe8eb27SJed Brown   } else {
27318be62a5SSatish Balay     ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
27418be62a5SSatish Balay   }
2758fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
2768fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
27718be62a5SSatish Balay   PetscFunctionReturn(0);
27818be62a5SSatish Balay }
27918be62a5SSatish Balay 
28018be62a5SSatish Balay #undef __FUNCT__
281ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell"
282f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
283ef66eb69SBarry Smith {
284ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
2858fe8eb27SJed Brown   PetscErrorCode ierr;
286b24ad042SBarry Smith 
287ef66eb69SBarry Smith   PetscFunctionBegin;
2888fe8eb27SJed Brown   if (shell->left || shell->right || shell->dshift) {
2898fe8eb27SJed Brown     if (!shell->dshift) {
2908fe8eb27SJed Brown       if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);}
2918fe8eb27SJed Brown       shell->dshift = shell->dshift_owned;
2928fe8eb27SJed Brown       ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr);
2938fe8eb27SJed Brown     } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
2948fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
2958fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
2968fe8eb27SJed Brown   } else shell->vshift += a;
2978fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
298ef66eb69SBarry Smith   PetscFunctionReturn(0);
299ef66eb69SBarry Smith }
300ef66eb69SBarry Smith 
301ef66eb69SBarry Smith #undef __FUNCT__
302ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell"
303f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
304ef66eb69SBarry Smith {
305ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3068fe8eb27SJed Brown   PetscErrorCode ierr;
307b24ad042SBarry Smith 
308ef66eb69SBarry Smith   PetscFunctionBegin;
309f4df32b1SMatthew Knepley   shell->vscale *= a;
3108fe8eb27SJed Brown   if (shell->dshift) {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
3118fe8eb27SJed Brown   else shell->vshift *= a;
3128fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
3138fe8eb27SJed Brown   PetscFunctionReturn(0);
31418be62a5SSatish Balay }
3158fe8eb27SJed Brown 
3168fe8eb27SJed Brown #undef __FUNCT__
3178fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell"
3188fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
3198fe8eb27SJed Brown {
3208fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)Y->data;
3218fe8eb27SJed Brown   PetscErrorCode ierr;
3228fe8eb27SJed Brown 
3238fe8eb27SJed Brown   PetscFunctionBegin;
3248fe8eb27SJed Brown   if (left) {
3258fe8eb27SJed Brown     if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);}
3268fe8eb27SJed Brown     if (shell->left) {ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);}
3278fe8eb27SJed Brown     else {
3288fe8eb27SJed Brown       shell->left = shell->left_owned;
3298fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
33018be62a5SSatish Balay     }
331ef66eb69SBarry Smith   }
3328fe8eb27SJed Brown   if (right) {
3338fe8eb27SJed Brown     if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);}
3348fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);}
3358fe8eb27SJed Brown     else {
3368fe8eb27SJed Brown       shell->right = shell->right_owned;
3378fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
3388fe8eb27SJed Brown     }
3398fe8eb27SJed Brown   }
3408fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
341ef66eb69SBarry Smith   PetscFunctionReturn(0);
342ef66eb69SBarry Smith }
343ef66eb69SBarry Smith 
344ef66eb69SBarry Smith #undef __FUNCT__
345ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell"
346dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
347ef66eb69SBarry Smith {
348ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
349ef66eb69SBarry Smith 
350ef66eb69SBarry Smith   PetscFunctionBegin;
3518fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
352ef66eb69SBarry Smith     shell->vshift      = 0.0;
353ef66eb69SBarry Smith     shell->vscale      = 1.0;
3548fe8eb27SJed Brown     shell->dshift      = PETSC_NULL;
3558fe8eb27SJed Brown     shell->left        = PETSC_NULL;
3568fe8eb27SJed Brown     shell->right       = PETSC_NULL;
3577fb7d96aSJed Brown     if (shell->mult) {
358ef66eb69SBarry Smith       Y->ops->mult = shell->mult;
3597fb7d96aSJed Brown       shell->mult  = PETSC_NULL;
3607fb7d96aSJed Brown     }
3617fb7d96aSJed Brown     if (shell->multtranspose) {
36218be62a5SSatish Balay       Y->ops->multtranspose = shell->multtranspose;
3637fb7d96aSJed Brown       shell->multtranspose  = PETSC_NULL;
3647fb7d96aSJed Brown     }
3657fb7d96aSJed Brown     if (shell->getdiagonal) {
36618be62a5SSatish Balay       Y->ops->getdiagonal = shell->getdiagonal;
3677fb7d96aSJed Brown       shell->getdiagonal  = PETSC_NULL;
3687fb7d96aSJed Brown     }
3697fb7d96aSJed Brown     shell->usingscaled = PETSC_FALSE;
370ef66eb69SBarry Smith   }
371ef66eb69SBarry Smith   PetscFunctionReturn(0);
372ef66eb69SBarry Smith }
373ef66eb69SBarry Smith 
37409573ac7SBarry Smith extern PetscErrorCode MatConvert_Shell(Mat, const MatType,MatReuse,Mat*);
375b951964fSBarry Smith 
37609dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
37720563c6bSBarry Smith        0,
37820563c6bSBarry Smith        0,
37920563c6bSBarry Smith        0,
38097304618SKris Buschelman /* 4*/ 0,
38120563c6bSBarry Smith        0,
382b951964fSBarry Smith        0,
383b951964fSBarry Smith        0,
384b951964fSBarry Smith        0,
385b951964fSBarry Smith        0,
38697304618SKris Buschelman /*10*/ 0,
387b951964fSBarry Smith        0,
388b951964fSBarry Smith        0,
389b951964fSBarry Smith        0,
390b951964fSBarry Smith        0,
39197304618SKris Buschelman /*15*/ 0,
392b951964fSBarry Smith        0,
393b951964fSBarry Smith        0,
3948fe8eb27SJed Brown        MatDiagonalScale_Shell,
395b951964fSBarry Smith        0,
39697304618SKris Buschelman /*20*/ 0,
397ef66eb69SBarry Smith        MatAssemblyEnd_Shell,
398b951964fSBarry Smith        0,
399b951964fSBarry Smith        0,
400d519adbfSMatthew Knepley /*24*/ 0,
401b951964fSBarry Smith        0,
402b951964fSBarry Smith        0,
403b951964fSBarry Smith        0,
404b951964fSBarry Smith        0,
405d519adbfSMatthew Knepley /*29*/ 0,
406b951964fSBarry Smith        0,
407273d9f13SBarry Smith        0,
408b951964fSBarry Smith        0,
409b951964fSBarry Smith        0,
410d519adbfSMatthew Knepley /*34*/ 0,
411b951964fSBarry Smith        0,
412b951964fSBarry Smith        0,
41309dc0095SBarry Smith        0,
41409dc0095SBarry Smith        0,
415d519adbfSMatthew Knepley /*39*/ 0,
41609dc0095SBarry Smith        0,
41709dc0095SBarry Smith        0,
41809dc0095SBarry Smith        0,
41909dc0095SBarry Smith        0,
420d519adbfSMatthew Knepley /*44*/ 0,
421ef66eb69SBarry Smith        MatScale_Shell,
422ef66eb69SBarry Smith        MatShift_Shell,
42309dc0095SBarry Smith        0,
42409dc0095SBarry Smith        0,
425f73d5cc4SBarry Smith /*49*/ 0,
42609dc0095SBarry Smith        0,
42709dc0095SBarry Smith        0,
42809dc0095SBarry Smith        0,
42909dc0095SBarry Smith        0,
430d519adbfSMatthew Knepley /*54*/ 0,
43109dc0095SBarry Smith        0,
43209dc0095SBarry Smith        0,
43309dc0095SBarry Smith        0,
43409dc0095SBarry Smith        0,
435d519adbfSMatthew Knepley /*59*/ 0,
436b9b97703SBarry Smith        MatDestroy_Shell,
43709dc0095SBarry Smith        0,
438357abbc8SBarry Smith        0,
439273d9f13SBarry Smith        0,
440d519adbfSMatthew Knepley /*64*/ 0,
441273d9f13SBarry Smith        0,
442273d9f13SBarry Smith        0,
443273d9f13SBarry Smith        0,
444273d9f13SBarry Smith        0,
445d519adbfSMatthew Knepley /*69*/ 0,
446273d9f13SBarry Smith        0,
447c87e5d42SMatthew Knepley        MatConvert_Shell,
448273d9f13SBarry Smith        0,
44997304618SKris Buschelman        0,
450d519adbfSMatthew Knepley /*74*/ 0,
45197304618SKris Buschelman        0,
45297304618SKris Buschelman        0,
45397304618SKris Buschelman        0,
45497304618SKris Buschelman        0,
455d519adbfSMatthew Knepley /*79*/ 0,
45697304618SKris Buschelman        0,
45797304618SKris Buschelman        0,
45897304618SKris Buschelman        0,
459865e5f61SKris Buschelman        0,
460d519adbfSMatthew Knepley /*84*/ 0,
461865e5f61SKris Buschelman        0,
462865e5f61SKris Buschelman        0,
463865e5f61SKris Buschelman        0,
464865e5f61SKris Buschelman        0,
465d519adbfSMatthew Knepley /*89*/ 0,
466865e5f61SKris Buschelman        0,
467865e5f61SKris Buschelman        0,
468865e5f61SKris Buschelman        0,
469865e5f61SKris Buschelman        0,
470d519adbfSMatthew Knepley /*94*/ 0,
471865e5f61SKris Buschelman        0,
472865e5f61SKris Buschelman        0,
473865e5f61SKris Buschelman        0};
474273d9f13SBarry Smith 
4750bad9183SKris Buschelman /*MC
476fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
4770bad9183SKris Buschelman 
4780bad9183SKris Buschelman   Level: advanced
4790bad9183SKris Buschelman 
4800bad9183SKris Buschelman .seealso: MatCreateShell
4810bad9183SKris Buschelman M*/
4820bad9183SKris Buschelman 
483273d9f13SBarry Smith EXTERN_C_BEGIN
484711e205bSSatish Balay #undef __FUNCT__
485711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell"
4867087cfbeSBarry Smith PetscErrorCode  MatCreate_Shell(Mat A)
487273d9f13SBarry Smith {
488273d9f13SBarry Smith   Mat_Shell      *b;
489dfbe8321SBarry Smith   PetscErrorCode ierr;
490273d9f13SBarry Smith 
491273d9f13SBarry Smith   PetscFunctionBegin;
492273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
493273d9f13SBarry Smith 
49438f2d2fdSLisandro Dalcin   ierr = PetscNewLog(A,Mat_Shell,&b);CHKERRQ(ierr);
495273d9f13SBarry Smith   A->data = (void*)b;
496273d9f13SBarry Smith 
49726283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
49826283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
499273d9f13SBarry Smith 
500273d9f13SBarry Smith   b->ctx           = 0;
501ef66eb69SBarry Smith   b->vshift        = 0.0;
502ef66eb69SBarry Smith   b->vscale        = 1.0;
503ef66eb69SBarry Smith   b->mult          = 0;
50418be62a5SSatish Balay   b->multtranspose = 0;
50518be62a5SSatish Balay   b->getdiagonal   = 0;
506273d9f13SBarry Smith   A->assembled     = PETSC_TRUE;
507210f0121SBarry Smith   A->preallocated  = PETSC_FALSE;
50817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
509273d9f13SBarry Smith   PetscFunctionReturn(0);
510273d9f13SBarry Smith }
511273d9f13SBarry Smith EXTERN_C_END
512e51e0e81SBarry Smith 
513711e205bSSatish Balay #undef __FUNCT__
514711e205bSSatish Balay #define __FUNCT__ "MatCreateShell"
5154b828684SBarry Smith /*@C
516052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
517ff756334SLois Curfman McInnes    private data storage format.
518e51e0e81SBarry Smith 
519c7fcc2eaSBarry Smith   Collective on MPI_Comm
520c7fcc2eaSBarry Smith 
521e51e0e81SBarry Smith    Input Parameters:
522c7fcc2eaSBarry Smith +  comm - MPI communicator
523c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
524c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
525c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
526c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
527c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
528e51e0e81SBarry Smith 
529ff756334SLois Curfman McInnes    Output Parameter:
53044cd7ae7SLois Curfman McInnes .  A - the matrix
531e51e0e81SBarry Smith 
532ff2fd236SBarry Smith    Level: advanced
533ff2fd236SBarry Smith 
534f39d1f56SLois Curfman McInnes   Usage:
5357b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
536f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
537c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
538f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
539f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
540f39d1f56SLois Curfman McInnes 
541ff756334SLois Curfman McInnes    Notes:
542ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
543ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
544ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
545e51e0e81SBarry Smith 
546f38a8d56SBarry Smith    Fortran Notes: The context can only be an integer or a PetscObject
547f38a8d56SBarry Smith       unfortunately it cannot be a Fortran array or derived type.
548f38a8d56SBarry Smith 
549f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
550f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
551645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
552645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
553645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
554645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
555645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
556645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
557645985a0SLois Curfman McInnes    For example,
558f39d1f56SLois Curfman McInnes 
559f39d1f56SLois Curfman McInnes $
560f39d1f56SLois Curfman McInnes $     Vec x, y
5617b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
562645985a0SLois Curfman McInnes $     Mat A
563f39d1f56SLois Curfman McInnes $
564c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
565c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
566f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
567c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
568c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
569c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
570645985a0SLois Curfman McInnes $     MatMult(A,x,y);
571645985a0SLois Curfman McInnes $     MatDestroy(A);
572f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
573645985a0SLois Curfman McInnes $
574e51e0e81SBarry Smith 
5750b627109SLois Curfman McInnes .keywords: matrix, shell, create
5760b627109SLois Curfman McInnes 
577ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
578e51e0e81SBarry Smith @*/
5797087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
580e51e0e81SBarry Smith {
581dfbe8321SBarry Smith   PetscErrorCode ierr;
582ed3cc1f0SBarry Smith 
5833a40ed3dSBarry Smith   PetscFunctionBegin;
584f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
585f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
586273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
587273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
5880fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
589273d9f13SBarry Smith   PetscFunctionReturn(0);
590c7fcc2eaSBarry Smith }
591c7fcc2eaSBarry Smith 
592711e205bSSatish Balay #undef __FUNCT__
593711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext"
594c6866cfdSSatish Balay /*@
595273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
596c7fcc2eaSBarry Smith 
5973f9fe445SBarry Smith    Logically Collective on Mat
598c7fcc2eaSBarry Smith 
599273d9f13SBarry Smith     Input Parameters:
600273d9f13SBarry Smith +   mat - the shell matrix
601273d9f13SBarry Smith -   ctx - the context
602273d9f13SBarry Smith 
603273d9f13SBarry Smith    Level: advanced
604273d9f13SBarry Smith 
605f38a8d56SBarry Smith    Fortran Notes: The context can only be an integer or a PetscObject
606f38a8d56SBarry Smith       unfortunately it cannot be a Fortran array or derived type.
607273d9f13SBarry Smith 
608273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
6090bc0a719SBarry Smith @*/
6107087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
611273d9f13SBarry Smith {
612273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
613dfbe8321SBarry Smith   PetscErrorCode ierr;
614ace3abfcSBarry Smith   PetscBool      flg;
615273d9f13SBarry Smith 
616273d9f13SBarry Smith   PetscFunctionBegin;
6170700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
618*251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
619273d9f13SBarry Smith   if (flg) {
620273d9f13SBarry Smith     shell->ctx = ctx;
621940183f3SBarry Smith   } else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
6223a40ed3dSBarry Smith   PetscFunctionReturn(0);
623e51e0e81SBarry Smith }
624e51e0e81SBarry Smith 
625711e205bSSatish Balay #undef __FUNCT__
626711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation"
627c16cb8f2SBarry Smith /*@C
6283a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
6293a3eedf2SBarry Smith                            a shell matrix.
630e51e0e81SBarry Smith 
6313f9fe445SBarry Smith    Logically Collective on Mat
632fee21e36SBarry Smith 
633c7fcc2eaSBarry Smith     Input Parameters:
634c7fcc2eaSBarry Smith +   mat - the shell matrix
635c7fcc2eaSBarry Smith .   op - the name of the operation
636c7fcc2eaSBarry Smith -   f - the function that provides the operation.
637c7fcc2eaSBarry Smith 
63815091d37SBarry Smith    Level: advanced
63915091d37SBarry Smith 
640fae171e0SBarry Smith     Usage:
641e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
642f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
643c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
6440b627109SLois Curfman McInnes 
645a62d957aSLois Curfman McInnes     Notes:
646e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
6471c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
648a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
6491c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
650a62d957aSLois Curfman McInnes 
65125296bd5SBarry Smith     All user-provided functions (execept for MATOP_DESTROY) should have the same calling
652deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
653deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
654deebb3c3SLois Curfman McInnes     routines, e.g.,
655a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
656a62d957aSLois Curfman McInnes 
6574aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
6584aa34b0aSBarry Smith     nonzero on failure.
6594aa34b0aSBarry Smith 
660a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
661a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
662a62d957aSLois Curfman McInnes     set by MatCreateShell().
663a62d957aSLois Curfman McInnes 
664501d9185SBarry Smith     Fortran Notes: For MatGetVecs() the user code should check if the input left or right matrix is -1 and in that case not
665501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
666501d9185SBarry Smith 
667a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
668a62d957aSLois Curfman McInnes 
669ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
670e51e0e81SBarry Smith @*/
6717087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
672e51e0e81SBarry Smith {
673dfbe8321SBarry Smith   PetscErrorCode ierr;
674ace3abfcSBarry Smith   PetscBool      flg;
675273d9f13SBarry Smith 
6763a40ed3dSBarry Smith   PetscFunctionBegin;
6770700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6781c1c02c0SLois Curfman McInnes   if (op == MATOP_DESTROY) {
679*251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
680273d9f13SBarry Smith     if (flg) {
681a62d957aSLois Curfman McInnes        Mat_Shell *shell = (Mat_Shell*)mat->data;
6826849ba73SBarry Smith        shell->destroy                 = (PetscErrorCode (*)(Mat)) f;
6836849ba73SBarry Smith     } else mat->ops->destroy          = (PetscErrorCode (*)(Mat)) f;
684a62d957aSLois Curfman McInnes   }
6856849ba73SBarry Smith   else if (op == MATOP_VIEW) mat->ops->view  = (PetscErrorCode (*)(Mat,PetscViewer)) f;
686c134de8dSSatish Balay   else                       (((void(**)(void))mat->ops)[op]) = f;
687a62d957aSLois Curfman McInnes 
6883a40ed3dSBarry Smith   PetscFunctionReturn(0);
689e51e0e81SBarry Smith }
690f0479e8cSBarry Smith 
691711e205bSSatish Balay #undef __FUNCT__
692711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation"
693d4bb536fSBarry Smith /*@C
694d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
695d4bb536fSBarry Smith 
696c7fcc2eaSBarry Smith     Not Collective
697c7fcc2eaSBarry Smith 
698d4bb536fSBarry Smith     Input Parameters:
699c7fcc2eaSBarry Smith +   mat - the shell matrix
700c7fcc2eaSBarry Smith -   op - the name of the operation
701d4bb536fSBarry Smith 
702d4bb536fSBarry Smith     Output Parameter:
703d4bb536fSBarry Smith .   f - the function that provides the operation.
704d4bb536fSBarry Smith 
70515091d37SBarry Smith     Level: advanced
70615091d37SBarry Smith 
707d4bb536fSBarry Smith     Notes:
708e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
709d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
710d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
711d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
712d4bb536fSBarry Smith 
713d4bb536fSBarry Smith     All user-provided functions have the same calling
714d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
715d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
716d4bb536fSBarry Smith     routines, e.g.,
717d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
718d4bb536fSBarry Smith 
719d4bb536fSBarry Smith     Within each user-defined routine, the user should call
720d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
721d4bb536fSBarry Smith     set by MatCreateShell().
722d4bb536fSBarry Smith 
723d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
724d4bb536fSBarry Smith 
725ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
726d4bb536fSBarry Smith @*/
7277087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
728d4bb536fSBarry Smith {
729dfbe8321SBarry Smith   PetscErrorCode ierr;
730ace3abfcSBarry Smith   PetscBool      flg;
731273d9f13SBarry Smith 
7323a40ed3dSBarry Smith   PetscFunctionBegin;
7330700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
734d4bb536fSBarry Smith   if (op == MATOP_DESTROY) {
735*251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
736273d9f13SBarry Smith     if (flg) {
737d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
738c134de8dSSatish Balay       *f = (void(*)(void))shell->destroy;
739c7fcc2eaSBarry Smith     } else {
740c134de8dSSatish Balay       *f = (void(*)(void))mat->ops->destroy;
741d4bb536fSBarry Smith     }
742c7fcc2eaSBarry Smith   } else if (op == MATOP_VIEW) {
743c134de8dSSatish Balay     *f = (void(*)(void))mat->ops->view;
744c7fcc2eaSBarry Smith   } else {
745c134de8dSSatish Balay     *f = (((void(**)(void))mat->ops)[op]);
746d4bb536fSBarry Smith   }
747d4bb536fSBarry Smith 
7483a40ed3dSBarry Smith   PetscFunctionReturn(0);
749d4bb536fSBarry Smith }
750d4bb536fSBarry Smith 
751