xref: /petsc/src/mat/impls/shell/shell.c (revision 7fb7d96aa659053597e5761399d97552ca6180df)
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 
8c6db04a5SJed Brown #include <private/matimpl.h>        /*I "petscmat.h" I*/
9c6db04a5SJed Brown #include <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;
898fe8eb27SJed Brown   if (!shell->left) {
908fe8eb27SJed Brown     *xx = x;
918fe8eb27SJed Brown   } else {
928fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
938fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
948fe8eb27SJed Brown     *xx = shell->left_work;
958fe8eb27SJed Brown   }
968fe8eb27SJed Brown   PetscFunctionReturn(0);
978fe8eb27SJed Brown }
988fe8eb27SJed Brown 
998fe8eb27SJed Brown #undef __FUNCT__
1008fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleRight"
1018fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
1028fe8eb27SJed Brown {
1038fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1048fe8eb27SJed Brown   PetscErrorCode ierr;
1058fe8eb27SJed Brown 
1068fe8eb27SJed Brown   PetscFunctionBegin;
1078fe8eb27SJed Brown   if (!shell->right) {
1088fe8eb27SJed Brown     *xx = x;
1098fe8eb27SJed Brown   } else {
1108fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
1118fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
1128fe8eb27SJed Brown     *xx = shell->right_work;
1138fe8eb27SJed Brown   }
1148fe8eb27SJed Brown   PetscFunctionReturn(0);
1158fe8eb27SJed Brown }
1168fe8eb27SJed Brown 
1178fe8eb27SJed Brown #undef __FUNCT__
1188fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleLeft"
1198fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
1208fe8eb27SJed Brown {
1218fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1228fe8eb27SJed Brown   PetscErrorCode ierr;
1238fe8eb27SJed Brown 
1248fe8eb27SJed Brown   PetscFunctionBegin;
1258fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
1268fe8eb27SJed Brown   PetscFunctionReturn(0);
1278fe8eb27SJed Brown }
1288fe8eb27SJed Brown 
1298fe8eb27SJed Brown #undef __FUNCT__
1308fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleRight"
1318fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
1328fe8eb27SJed Brown {
1338fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1348fe8eb27SJed Brown   PetscErrorCode ierr;
1358fe8eb27SJed Brown 
1368fe8eb27SJed Brown   PetscFunctionBegin;
1378fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
1388fe8eb27SJed Brown   PetscFunctionReturn(0);
1398fe8eb27SJed Brown }
1408fe8eb27SJed Brown 
1418fe8eb27SJed Brown #undef __FUNCT__
1428fe8eb27SJed Brown #define __FUNCT__ "MatShellShiftAndScale"
1438fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
1448fe8eb27SJed Brown {
1458fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)A->data;
1468fe8eb27SJed Brown   PetscErrorCode ierr;
1478fe8eb27SJed Brown 
1488fe8eb27SJed Brown   PetscFunctionBegin;
1498fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
1508fe8eb27SJed Brown     PetscInt i,m;
1518fe8eb27SJed Brown     const PetscScalar *x,*d;
1528fe8eb27SJed Brown     PetscScalar *y;
1538fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
1548fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1558fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1568fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1578fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
1588fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1598fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1608fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
1618fe8eb27SJed Brown   } else {
1628fe8eb27SJed Brown     ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr);
1638fe8eb27SJed Brown   }
1648fe8eb27SJed Brown   PetscFunctionReturn(0);
1658fe8eb27SJed Brown }
166e51e0e81SBarry Smith 
167711e205bSSatish Balay #undef __FUNCT__
168711e205bSSatish Balay #define __FUNCT__ "MatShellGetContext"
16935d520bfSBarry Smith /*@C
170a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
171b4fd4287SBarry Smith 
172c7fcc2eaSBarry Smith     Not Collective
173c7fcc2eaSBarry Smith 
174b4fd4287SBarry Smith     Input Parameter:
175b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
176b4fd4287SBarry Smith 
177b4fd4287SBarry Smith     Output Parameter:
178b4fd4287SBarry Smith .   ctx - the user provided context
179b4fd4287SBarry Smith 
18015091d37SBarry Smith     Level: advanced
18115091d37SBarry Smith 
182a62d957aSLois Curfman McInnes     Notes:
183a62d957aSLois Curfman McInnes     This routine is intended for use within various shell matrix routines,
184a62d957aSLois Curfman McInnes     as set with MatShellSetOperation().
185a62d957aSLois Curfman McInnes 
186a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context
187a62d957aSLois Curfman McInnes 
188ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
1890bc0a719SBarry Smith @*/
1908fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
191b4fd4287SBarry Smith {
192dfbe8321SBarry Smith   PetscErrorCode ierr;
193ace3abfcSBarry Smith   PetscBool      flg;
194273d9f13SBarry Smith 
1953a40ed3dSBarry Smith   PetscFunctionBegin;
1960700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1974482741eSBarry Smith   PetscValidPointer(ctx,2);
198273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1998fe8eb27SJed Brown   if (!flg) *(void**)ctx = 0;
2008fe8eb27SJed Brown   else      *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
2013a40ed3dSBarry Smith   PetscFunctionReturn(0);
202b4fd4287SBarry Smith }
203b4fd4287SBarry Smith 
204711e205bSSatish Balay #undef __FUNCT__
205711e205bSSatish Balay #define __FUNCT__ "MatDestroy_Shell"
206dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
207e51e0e81SBarry Smith {
208dfbe8321SBarry Smith   PetscErrorCode ierr;
209bf0cc555SLisandro Dalcin   Mat_Shell      *shell = (Mat_Shell*)mat->data;
210ed3cc1f0SBarry Smith 
2113a40ed3dSBarry Smith   PetscFunctionBegin;
212bf0cc555SLisandro Dalcin   if (shell->destroy) {
213bf0cc555SLisandro Dalcin     ierr = (*shell->destroy)(mat);CHKERRQ(ierr);
214bf0cc555SLisandro Dalcin   }
2158fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr);
2168fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr);
2178fe8eb27SJed Brown   ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr);
2188fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
2198fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
220bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
2213a40ed3dSBarry Smith   PetscFunctionReturn(0);
222e51e0e81SBarry Smith }
223e51e0e81SBarry Smith 
224711e205bSSatish Balay #undef __FUNCT__
225ef66eb69SBarry Smith #define __FUNCT__ "MatMult_Shell"
226dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
227ef66eb69SBarry Smith {
228ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)A->data;
229dfbe8321SBarry Smith   PetscErrorCode ierr;
2308fe8eb27SJed Brown   Vec            xx;
231ef66eb69SBarry Smith 
232ef66eb69SBarry Smith   PetscFunctionBegin;
2338fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
2348fe8eb27SJed Brown   ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr);
2358fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2368fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
237ef66eb69SBarry Smith   PetscFunctionReturn(0);
238ef66eb69SBarry Smith }
239ef66eb69SBarry Smith 
240ef66eb69SBarry Smith #undef __FUNCT__
24118be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell"
24218be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
24318be62a5SSatish Balay {
24418be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
24518be62a5SSatish Balay   PetscErrorCode ierr;
2468fe8eb27SJed Brown   Vec            xx;
24718be62a5SSatish Balay 
24818be62a5SSatish Balay   PetscFunctionBegin;
2498fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
2508fe8eb27SJed Brown   ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr);
2518fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2528fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
25318be62a5SSatish Balay   PetscFunctionReturn(0);
25418be62a5SSatish Balay }
25518be62a5SSatish Balay 
25618be62a5SSatish Balay #undef __FUNCT__
25718be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell"
25818be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
25918be62a5SSatish Balay {
26018be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
26118be62a5SSatish Balay   PetscErrorCode ierr;
26218be62a5SSatish Balay 
26318be62a5SSatish Balay   PetscFunctionBegin;
26418be62a5SSatish Balay   ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr);
26518be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
2668fe8eb27SJed Brown   if (shell->dshift) {
2678fe8eb27SJed Brown     ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr);
2688fe8eb27SJed Brown   } else {
26918be62a5SSatish Balay     ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
27018be62a5SSatish Balay   }
2718fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
2728fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
27318be62a5SSatish Balay   PetscFunctionReturn(0);
27418be62a5SSatish Balay }
27518be62a5SSatish Balay 
27618be62a5SSatish Balay #undef __FUNCT__
277ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell"
278f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
279ef66eb69SBarry Smith {
280ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
2818fe8eb27SJed Brown   PetscErrorCode ierr;
282b24ad042SBarry Smith 
283ef66eb69SBarry Smith   PetscFunctionBegin;
2848fe8eb27SJed Brown   if (shell->left || shell->right || shell->dshift) {
2858fe8eb27SJed Brown     if (!shell->dshift) {
2868fe8eb27SJed Brown       if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);}
2878fe8eb27SJed Brown       shell->dshift = shell->dshift_owned;
2888fe8eb27SJed Brown       ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr);
2898fe8eb27SJed Brown     } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
2908fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
2918fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
2928fe8eb27SJed Brown   } else shell->vshift += a;
2938fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
294ef66eb69SBarry Smith   PetscFunctionReturn(0);
295ef66eb69SBarry Smith }
296ef66eb69SBarry Smith 
297ef66eb69SBarry Smith #undef __FUNCT__
298ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell"
299f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
300ef66eb69SBarry Smith {
301ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3028fe8eb27SJed Brown   PetscErrorCode ierr;
303b24ad042SBarry Smith 
304ef66eb69SBarry Smith   PetscFunctionBegin;
305f4df32b1SMatthew Knepley   shell->vscale *= a;
3068fe8eb27SJed Brown   if (shell->dshift) {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
3078fe8eb27SJed Brown   else shell->vshift *= a;
3088fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
3098fe8eb27SJed Brown   PetscFunctionReturn(0);
31018be62a5SSatish Balay }
3118fe8eb27SJed Brown 
3128fe8eb27SJed Brown #undef __FUNCT__
3138fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell"
3148fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
3158fe8eb27SJed Brown {
3168fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)Y->data;
3178fe8eb27SJed Brown   PetscErrorCode ierr;
3188fe8eb27SJed Brown 
3198fe8eb27SJed Brown   PetscFunctionBegin;
3208fe8eb27SJed Brown   if (left) {
3218fe8eb27SJed Brown     if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);}
3228fe8eb27SJed Brown     if (shell->left) {ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);}
3238fe8eb27SJed Brown     else {
3248fe8eb27SJed Brown       shell->left = shell->left_owned;
3258fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
32618be62a5SSatish Balay     }
327ef66eb69SBarry Smith   }
3288fe8eb27SJed Brown   if (right) {
3298fe8eb27SJed Brown     if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);}
3308fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);}
3318fe8eb27SJed Brown     else {
3328fe8eb27SJed Brown       shell->right = shell->right_owned;
3338fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
3348fe8eb27SJed Brown     }
3358fe8eb27SJed Brown   }
3368fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
337ef66eb69SBarry Smith   PetscFunctionReturn(0);
338ef66eb69SBarry Smith }
339ef66eb69SBarry Smith 
340ef66eb69SBarry Smith #undef __FUNCT__
341ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell"
342dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
343ef66eb69SBarry Smith {
344ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
345ef66eb69SBarry Smith 
346ef66eb69SBarry Smith   PetscFunctionBegin;
3478fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
348ef66eb69SBarry Smith     shell->vshift      = 0.0;
349ef66eb69SBarry Smith     shell->vscale      = 1.0;
3508fe8eb27SJed Brown     shell->dshift      = PETSC_NULL;
3518fe8eb27SJed Brown     shell->left        = PETSC_NULL;
3528fe8eb27SJed Brown     shell->right       = PETSC_NULL;
353*7fb7d96aSJed Brown     if (shell->mult) {
354ef66eb69SBarry Smith       Y->ops->mult = shell->mult;
355*7fb7d96aSJed Brown       shell->mult  = PETSC_NULL;
356*7fb7d96aSJed Brown     }
357*7fb7d96aSJed Brown     if (shell->multtranspose) {
35818be62a5SSatish Balay       Y->ops->multtranspose = shell->multtranspose;
359*7fb7d96aSJed Brown       shell->multtranspose  = PETSC_NULL;
360*7fb7d96aSJed Brown     }
361*7fb7d96aSJed Brown     if (shell->getdiagonal) {
36218be62a5SSatish Balay       Y->ops->getdiagonal = shell->getdiagonal;
363*7fb7d96aSJed Brown       shell->getdiagonal  = PETSC_NULL;
364*7fb7d96aSJed Brown     }
365*7fb7d96aSJed Brown     shell->usingscaled = PETSC_FALSE;
366ef66eb69SBarry Smith   }
367ef66eb69SBarry Smith   PetscFunctionReturn(0);
368ef66eb69SBarry Smith }
369ef66eb69SBarry Smith 
37009573ac7SBarry Smith extern PetscErrorCode MatConvert_Shell(Mat, const MatType,MatReuse,Mat*);
371b951964fSBarry Smith 
372521d7252SBarry Smith #undef __FUNCT__
373521d7252SBarry Smith #define __FUNCT__ "MatSetBlockSize_Shell"
374521d7252SBarry Smith PetscErrorCode MatSetBlockSize_Shell(Mat A,PetscInt bs)
375521d7252SBarry Smith {
37641c166b1SJed Brown   PetscErrorCode ierr;
37741c166b1SJed Brown 
378521d7252SBarry Smith   PetscFunctionBegin;
37941c166b1SJed Brown   ierr = PetscLayoutSetBlockSize(A->rmap,bs);CHKERRQ(ierr);
38041c166b1SJed Brown   ierr = PetscLayoutSetBlockSize(A->cmap,bs);CHKERRQ(ierr);
381521d7252SBarry Smith   PetscFunctionReturn(0);
382521d7252SBarry Smith }
383521d7252SBarry Smith 
38409dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
38520563c6bSBarry Smith        0,
38620563c6bSBarry Smith        0,
38720563c6bSBarry Smith        0,
38897304618SKris Buschelman /* 4*/ 0,
38920563c6bSBarry Smith        0,
390b951964fSBarry Smith        0,
391b951964fSBarry Smith        0,
392b951964fSBarry Smith        0,
393b951964fSBarry Smith        0,
39497304618SKris Buschelman /*10*/ 0,
395b951964fSBarry Smith        0,
396b951964fSBarry Smith        0,
397b951964fSBarry Smith        0,
398b951964fSBarry Smith        0,
39997304618SKris Buschelman /*15*/ 0,
400b951964fSBarry Smith        0,
401b951964fSBarry Smith        0,
4028fe8eb27SJed Brown        MatDiagonalScale_Shell,
403b951964fSBarry Smith        0,
40497304618SKris Buschelman /*20*/ 0,
405ef66eb69SBarry Smith        MatAssemblyEnd_Shell,
406b951964fSBarry Smith        0,
407b951964fSBarry Smith        0,
408d519adbfSMatthew Knepley /*24*/ 0,
409b951964fSBarry Smith        0,
410b951964fSBarry Smith        0,
411b951964fSBarry Smith        0,
412b951964fSBarry Smith        0,
413d519adbfSMatthew Knepley /*29*/ 0,
414b951964fSBarry Smith        0,
415273d9f13SBarry Smith        0,
416b951964fSBarry Smith        0,
417b951964fSBarry Smith        0,
418d519adbfSMatthew Knepley /*34*/ 0,
419b951964fSBarry Smith        0,
420b951964fSBarry Smith        0,
42109dc0095SBarry Smith        0,
42209dc0095SBarry Smith        0,
423d519adbfSMatthew Knepley /*39*/ 0,
42409dc0095SBarry Smith        0,
42509dc0095SBarry Smith        0,
42609dc0095SBarry Smith        0,
42709dc0095SBarry Smith        0,
428d519adbfSMatthew Knepley /*44*/ 0,
429ef66eb69SBarry Smith        MatScale_Shell,
430ef66eb69SBarry Smith        MatShift_Shell,
43109dc0095SBarry Smith        0,
43209dc0095SBarry Smith        0,
433d519adbfSMatthew Knepley /*49*/ MatSetBlockSize_Shell,
43409dc0095SBarry Smith        0,
43509dc0095SBarry Smith        0,
43609dc0095SBarry Smith        0,
43709dc0095SBarry Smith        0,
438d519adbfSMatthew Knepley /*54*/ 0,
43909dc0095SBarry Smith        0,
44009dc0095SBarry Smith        0,
44109dc0095SBarry Smith        0,
44209dc0095SBarry Smith        0,
443d519adbfSMatthew Knepley /*59*/ 0,
444b9b97703SBarry Smith        MatDestroy_Shell,
44509dc0095SBarry Smith        0,
446357abbc8SBarry Smith        0,
447273d9f13SBarry Smith        0,
448d519adbfSMatthew Knepley /*64*/ 0,
449273d9f13SBarry Smith        0,
450273d9f13SBarry Smith        0,
451273d9f13SBarry Smith        0,
452273d9f13SBarry Smith        0,
453d519adbfSMatthew Knepley /*69*/ 0,
454273d9f13SBarry Smith        0,
455c87e5d42SMatthew Knepley        MatConvert_Shell,
456273d9f13SBarry Smith        0,
45797304618SKris Buschelman        0,
458d519adbfSMatthew Knepley /*74*/ 0,
45997304618SKris Buschelman        0,
46097304618SKris Buschelman        0,
46197304618SKris Buschelman        0,
46297304618SKris Buschelman        0,
463d519adbfSMatthew Knepley /*79*/ 0,
46497304618SKris Buschelman        0,
46597304618SKris Buschelman        0,
46697304618SKris Buschelman        0,
467865e5f61SKris Buschelman        0,
468d519adbfSMatthew Knepley /*84*/ 0,
469865e5f61SKris Buschelman        0,
470865e5f61SKris Buschelman        0,
471865e5f61SKris Buschelman        0,
472865e5f61SKris Buschelman        0,
473d519adbfSMatthew Knepley /*89*/ 0,
474865e5f61SKris Buschelman        0,
475865e5f61SKris Buschelman        0,
476865e5f61SKris Buschelman        0,
477865e5f61SKris Buschelman        0,
478d519adbfSMatthew Knepley /*94*/ 0,
479865e5f61SKris Buschelman        0,
480865e5f61SKris Buschelman        0,
481865e5f61SKris Buschelman        0};
482273d9f13SBarry Smith 
4830bad9183SKris Buschelman /*MC
484fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
4850bad9183SKris Buschelman 
4860bad9183SKris Buschelman   Level: advanced
4870bad9183SKris Buschelman 
4880bad9183SKris Buschelman .seealso: MatCreateShell
4890bad9183SKris Buschelman M*/
4900bad9183SKris Buschelman 
491273d9f13SBarry Smith EXTERN_C_BEGIN
492711e205bSSatish Balay #undef __FUNCT__
493711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell"
4947087cfbeSBarry Smith PetscErrorCode  MatCreate_Shell(Mat A)
495273d9f13SBarry Smith {
496273d9f13SBarry Smith   Mat_Shell      *b;
497dfbe8321SBarry Smith   PetscErrorCode ierr;
498273d9f13SBarry Smith 
499273d9f13SBarry Smith   PetscFunctionBegin;
500273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
501273d9f13SBarry Smith 
50238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(A,Mat_Shell,&b);CHKERRQ(ierr);
503273d9f13SBarry Smith   A->data = (void*)b;
504273d9f13SBarry Smith 
50526283091SBarry Smith   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
50626283091SBarry Smith   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
50726283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
50826283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
509273d9f13SBarry Smith 
510273d9f13SBarry Smith   b->ctx           = 0;
511ef66eb69SBarry Smith   b->vshift        = 0.0;
512ef66eb69SBarry Smith   b->vscale        = 1.0;
513ef66eb69SBarry Smith   b->mult          = 0;
51418be62a5SSatish Balay   b->multtranspose = 0;
51518be62a5SSatish Balay   b->getdiagonal   = 0;
516273d9f13SBarry Smith   A->assembled     = PETSC_TRUE;
517210f0121SBarry Smith   A->preallocated  = PETSC_FALSE;
51817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
519273d9f13SBarry Smith   PetscFunctionReturn(0);
520273d9f13SBarry Smith }
521273d9f13SBarry Smith EXTERN_C_END
522e51e0e81SBarry Smith 
523711e205bSSatish Balay #undef __FUNCT__
524711e205bSSatish Balay #define __FUNCT__ "MatCreateShell"
5254b828684SBarry Smith /*@C
526052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
527ff756334SLois Curfman McInnes    private data storage format.
528e51e0e81SBarry Smith 
529c7fcc2eaSBarry Smith   Collective on MPI_Comm
530c7fcc2eaSBarry Smith 
531e51e0e81SBarry Smith    Input Parameters:
532c7fcc2eaSBarry Smith +  comm - MPI communicator
533c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
534c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
535c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
536c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
537c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
538e51e0e81SBarry Smith 
539ff756334SLois Curfman McInnes    Output Parameter:
54044cd7ae7SLois Curfman McInnes .  A - the matrix
541e51e0e81SBarry Smith 
542ff2fd236SBarry Smith    Level: advanced
543ff2fd236SBarry Smith 
544f39d1f56SLois Curfman McInnes   Usage:
5457b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
546f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
547c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
548f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
549f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
550f39d1f56SLois Curfman McInnes 
551ff756334SLois Curfman McInnes    Notes:
552ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
553ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
554ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
555e51e0e81SBarry Smith 
556f38a8d56SBarry Smith    Fortran Notes: The context can only be an integer or a PetscObject
557f38a8d56SBarry Smith       unfortunately it cannot be a Fortran array or derived type.
558f38a8d56SBarry Smith 
559f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
560f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
561645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
562645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
563645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
564645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
565645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
566645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
567645985a0SLois Curfman McInnes    For example,
568f39d1f56SLois Curfman McInnes 
569f39d1f56SLois Curfman McInnes $
570f39d1f56SLois Curfman McInnes $     Vec x, y
5717b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
572645985a0SLois Curfman McInnes $     Mat A
573f39d1f56SLois Curfman McInnes $
574c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
575c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
576f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
577c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
578c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
579c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
580645985a0SLois Curfman McInnes $     MatMult(A,x,y);
581645985a0SLois Curfman McInnes $     MatDestroy(A);
582f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
583645985a0SLois Curfman McInnes $
584e51e0e81SBarry Smith 
5850b627109SLois Curfman McInnes .keywords: matrix, shell, create
5860b627109SLois Curfman McInnes 
587ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
588e51e0e81SBarry Smith @*/
5897087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
590e51e0e81SBarry Smith {
591dfbe8321SBarry Smith   PetscErrorCode ierr;
592ed3cc1f0SBarry Smith 
5933a40ed3dSBarry Smith   PetscFunctionBegin;
594f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
595f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
596899cda47SBarry Smith 
597273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
598273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
599273d9f13SBarry Smith   PetscFunctionReturn(0);
600c7fcc2eaSBarry Smith }
601c7fcc2eaSBarry Smith 
602711e205bSSatish Balay #undef __FUNCT__
603711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext"
604c6866cfdSSatish Balay /*@
605273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
606c7fcc2eaSBarry Smith 
6073f9fe445SBarry Smith    Logically Collective on Mat
608c7fcc2eaSBarry Smith 
609273d9f13SBarry Smith     Input Parameters:
610273d9f13SBarry Smith +   mat - the shell matrix
611273d9f13SBarry Smith -   ctx - the context
612273d9f13SBarry Smith 
613273d9f13SBarry Smith    Level: advanced
614273d9f13SBarry Smith 
615f38a8d56SBarry Smith    Fortran Notes: The context can only be an integer or a PetscObject
616f38a8d56SBarry Smith       unfortunately it cannot be a Fortran array or derived type.
617273d9f13SBarry Smith 
618273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
6190bc0a719SBarry Smith @*/
6207087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
621273d9f13SBarry Smith {
622273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
623dfbe8321SBarry Smith   PetscErrorCode ierr;
624ace3abfcSBarry Smith   PetscBool      flg;
625273d9f13SBarry Smith 
626273d9f13SBarry Smith   PetscFunctionBegin;
6270700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
628273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
629273d9f13SBarry Smith   if (flg) {
630273d9f13SBarry Smith     shell->ctx = ctx;
631273d9f13SBarry Smith   }
6323a40ed3dSBarry Smith   PetscFunctionReturn(0);
633e51e0e81SBarry Smith }
634e51e0e81SBarry Smith 
635711e205bSSatish Balay #undef __FUNCT__
636711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation"
637c16cb8f2SBarry Smith /*@C
6383a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
6393a3eedf2SBarry Smith                            a shell matrix.
640e51e0e81SBarry Smith 
6413f9fe445SBarry Smith    Logically Collective on Mat
642fee21e36SBarry Smith 
643c7fcc2eaSBarry Smith     Input Parameters:
644c7fcc2eaSBarry Smith +   mat - the shell matrix
645c7fcc2eaSBarry Smith .   op - the name of the operation
646c7fcc2eaSBarry Smith -   f - the function that provides the operation.
647c7fcc2eaSBarry Smith 
64815091d37SBarry Smith    Level: advanced
64915091d37SBarry Smith 
650fae171e0SBarry Smith     Usage:
651e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
652f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
653c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
6540b627109SLois Curfman McInnes 
655a62d957aSLois Curfman McInnes     Notes:
656e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
6571c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
658a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
6591c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
660a62d957aSLois Curfman McInnes 
661a62d957aSLois Curfman McInnes     All user-provided functions should have the same calling
662deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
663deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
664deebb3c3SLois Curfman McInnes     routines, e.g.,
665a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
666a62d957aSLois Curfman McInnes 
6674aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
6684aa34b0aSBarry Smith     nonzero on failure.
6694aa34b0aSBarry Smith 
670a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
671a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
672a62d957aSLois Curfman McInnes     set by MatCreateShell().
673a62d957aSLois Curfman McInnes 
674501d9185SBarry Smith     Fortran Notes: For MatGetVecs() the user code should check if the input left or right matrix is -1 and in that case not
675501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
676501d9185SBarry Smith 
677a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
678a62d957aSLois Curfman McInnes 
679ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
680e51e0e81SBarry Smith @*/
6817087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
682e51e0e81SBarry Smith {
683dfbe8321SBarry Smith   PetscErrorCode ierr;
684ace3abfcSBarry Smith   PetscBool      flg;
685273d9f13SBarry Smith 
6863a40ed3dSBarry Smith   PetscFunctionBegin;
6870700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6881c1c02c0SLois Curfman McInnes   if (op == MATOP_DESTROY) {
689273d9f13SBarry Smith     ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
690273d9f13SBarry Smith     if (flg) {
691a62d957aSLois Curfman McInnes        Mat_Shell *shell = (Mat_Shell*)mat->data;
6926849ba73SBarry Smith        shell->destroy                 = (PetscErrorCode (*)(Mat)) f;
6936849ba73SBarry Smith     } else mat->ops->destroy          = (PetscErrorCode (*)(Mat)) f;
694a62d957aSLois Curfman McInnes   }
6956849ba73SBarry Smith   else if (op == MATOP_VIEW) mat->ops->view  = (PetscErrorCode (*)(Mat,PetscViewer)) f;
696c134de8dSSatish Balay   else                       (((void(**)(void))mat->ops)[op]) = f;
697a62d957aSLois Curfman McInnes 
6983a40ed3dSBarry Smith   PetscFunctionReturn(0);
699e51e0e81SBarry Smith }
700f0479e8cSBarry Smith 
701711e205bSSatish Balay #undef __FUNCT__
702711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation"
703d4bb536fSBarry Smith /*@C
704d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
705d4bb536fSBarry Smith 
706c7fcc2eaSBarry Smith     Not Collective
707c7fcc2eaSBarry Smith 
708d4bb536fSBarry Smith     Input Parameters:
709c7fcc2eaSBarry Smith +   mat - the shell matrix
710c7fcc2eaSBarry Smith -   op - the name of the operation
711d4bb536fSBarry Smith 
712d4bb536fSBarry Smith     Output Parameter:
713d4bb536fSBarry Smith .   f - the function that provides the operation.
714d4bb536fSBarry Smith 
71515091d37SBarry Smith     Level: advanced
71615091d37SBarry Smith 
717d4bb536fSBarry Smith     Notes:
718e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
719d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
720d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
721d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
722d4bb536fSBarry Smith 
723d4bb536fSBarry Smith     All user-provided functions have the same calling
724d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
725d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
726d4bb536fSBarry Smith     routines, e.g.,
727d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
728d4bb536fSBarry Smith 
729d4bb536fSBarry Smith     Within each user-defined routine, the user should call
730d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
731d4bb536fSBarry Smith     set by MatCreateShell().
732d4bb536fSBarry Smith 
733d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
734d4bb536fSBarry Smith 
735ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
736d4bb536fSBarry Smith @*/
7377087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
738d4bb536fSBarry Smith {
739dfbe8321SBarry Smith   PetscErrorCode ierr;
740ace3abfcSBarry Smith   PetscBool      flg;
741273d9f13SBarry Smith 
7423a40ed3dSBarry Smith   PetscFunctionBegin;
7430700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
744d4bb536fSBarry Smith   if (op == MATOP_DESTROY) {
745273d9f13SBarry Smith     ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
746273d9f13SBarry Smith     if (flg) {
747d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
748c134de8dSSatish Balay       *f = (void(*)(void))shell->destroy;
749c7fcc2eaSBarry Smith     } else {
750c134de8dSSatish Balay       *f = (void(*)(void))mat->ops->destroy;
751d4bb536fSBarry Smith     }
752c7fcc2eaSBarry Smith   } else if (op == MATOP_VIEW) {
753c134de8dSSatish Balay     *f = (void(*)(void))mat->ops->view;
754c7fcc2eaSBarry Smith   } else {
755c134de8dSSatish Balay     *f = (((void(**)(void))mat->ops)[op]);
756d4bb536fSBarry Smith   }
757d4bb536fSBarry Smith 
7583a40ed3dSBarry Smith   PetscFunctionReturn(0);
759d4bb536fSBarry Smith }
760d4bb536fSBarry Smith 
761