xref: /petsc/src/mat/impls/shell/shell.c (revision e3f487b05763f6cc9a19c71bd33970cdca741ac0)
1be1d678aSKris Buschelman 
2e51e0e81SBarry Smith /*
320563c6bSBarry Smith    This provides a simple shell for Fortran (and C programmers) to
420563c6bSBarry Smith   create a very simple matrix class for use with KSP without coding
5ed3cc1f0SBarry Smith   much of anything.
6e51e0e81SBarry Smith */
7e51e0e81SBarry Smith 
8af0996ceSBarry Smith #include <petsc/private/matimpl.h>        /*I "petscmat.h" I*/
9af0996ceSBarry Smith #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);
162205254eSKarl Rupp 
17ef66eb69SBarry Smith   PetscScalar vscale,vshift;
188fe8eb27SJed Brown   Vec         dshift;
198fe8eb27SJed Brown   Vec         left,right;
208fe8eb27SJed Brown   Vec         dshift_owned,left_owned,right_owned;
218fe8eb27SJed Brown   Vec         left_work,right_work;
225edf6516SJed Brown   Vec         left_add_work,right_add_work;
238fe8eb27SJed Brown   PetscBool   usingscaled;
2420563c6bSBarry Smith   void        *ctx;
2588cf3e7dSBarry Smith } Mat_Shell;
268fe8eb27SJed Brown /*
278fe8eb27SJed Brown  The most general expression for the matrix is
288fe8eb27SJed Brown 
298fe8eb27SJed Brown  S = L (a A + B) R
308fe8eb27SJed Brown 
318fe8eb27SJed Brown  where
328fe8eb27SJed Brown  A is the matrix defined by the user's function
338fe8eb27SJed Brown  a is a scalar multiple
348fe8eb27SJed Brown  L is left scaling
358fe8eb27SJed Brown  R is right scaling
368fe8eb27SJed Brown  B is a diagonal shift defined by
378fe8eb27SJed Brown    diag(dshift) if the vector dshift is non-NULL
388fe8eb27SJed Brown    vscale*identity otherwise
398fe8eb27SJed Brown 
408fe8eb27SJed Brown  The following identities apply:
418fe8eb27SJed Brown 
428fe8eb27SJed Brown  Scale by c:
438fe8eb27SJed Brown   c [L (a A + B) R] = L [(a c) A + c B] R
448fe8eb27SJed Brown 
458fe8eb27SJed Brown  Shift by c:
468fe8eb27SJed Brown   [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R
478fe8eb27SJed Brown 
488fe8eb27SJed Brown  Diagonal scaling is achieved by simply multiplying with existing L and R vectors
498fe8eb27SJed Brown 
508fe8eb27SJed Brown  In the data structure:
518fe8eb27SJed Brown 
528fe8eb27SJed Brown  vscale=1.0  means no special scaling will be applied
538fe8eb27SJed Brown  dshift=NULL means a constant diagonal shift (fall back to vshift)
548fe8eb27SJed Brown  vshift=0.0  means no constant diagonal shift, note that vshift is only used if dshift is NULL
558fe8eb27SJed Brown */
568fe8eb27SJed Brown 
578fe8eb27SJed Brown static PetscErrorCode MatMult_Shell(Mat,Vec,Vec);
588fe8eb27SJed Brown static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec);
598fe8eb27SJed Brown static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec);
608fe8eb27SJed Brown 
618fe8eb27SJed Brown #undef __FUNCT__
628fe8eb27SJed Brown #define __FUNCT__ "MatShellUseScaledMethods"
638fe8eb27SJed Brown static PetscErrorCode MatShellUseScaledMethods(Mat Y)
648fe8eb27SJed Brown {
658fe8eb27SJed Brown   Mat_Shell *shell = (Mat_Shell*)Y->data;
668fe8eb27SJed Brown 
678fe8eb27SJed Brown   PetscFunctionBegin;
688fe8eb27SJed Brown   if (shell->usingscaled) PetscFunctionReturn(0);
698fe8eb27SJed Brown   shell->mult  = Y->ops->mult;
708fe8eb27SJed Brown   Y->ops->mult = MatMult_Shell;
718fe8eb27SJed Brown   if (Y->ops->multtranspose) {
728fe8eb27SJed Brown     shell->multtranspose  = Y->ops->multtranspose;
738fe8eb27SJed Brown     Y->ops->multtranspose = MatMultTranspose_Shell;
748fe8eb27SJed Brown   }
758fe8eb27SJed Brown   if (Y->ops->getdiagonal) {
768fe8eb27SJed Brown     shell->getdiagonal  = Y->ops->getdiagonal;
778fe8eb27SJed Brown     Y->ops->getdiagonal = MatGetDiagonal_Shell;
788fe8eb27SJed Brown   }
798fe8eb27SJed Brown   shell->usingscaled = PETSC_TRUE;
808fe8eb27SJed Brown   PetscFunctionReturn(0);
818fe8eb27SJed Brown }
828fe8eb27SJed Brown 
838fe8eb27SJed Brown #undef __FUNCT__
848fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleLeft"
858fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
868fe8eb27SJed Brown {
878fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
888fe8eb27SJed Brown   PetscErrorCode ierr;
898fe8eb27SJed Brown 
908fe8eb27SJed Brown   PetscFunctionBegin;
910298fd71SBarry Smith   *xx = NULL;
928fe8eb27SJed Brown   if (!shell->left) {
938fe8eb27SJed Brown     *xx = x;
948fe8eb27SJed Brown   } else {
958fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
968fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
978fe8eb27SJed Brown     *xx  = shell->left_work;
988fe8eb27SJed Brown   }
998fe8eb27SJed Brown   PetscFunctionReturn(0);
1008fe8eb27SJed Brown }
1018fe8eb27SJed Brown 
1028fe8eb27SJed Brown #undef __FUNCT__
1038fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleRight"
1048fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
1058fe8eb27SJed Brown {
1068fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1078fe8eb27SJed Brown   PetscErrorCode ierr;
1088fe8eb27SJed Brown 
1098fe8eb27SJed Brown   PetscFunctionBegin;
1100298fd71SBarry Smith   *xx = NULL;
1118fe8eb27SJed Brown   if (!shell->right) {
1128fe8eb27SJed Brown     *xx = x;
1138fe8eb27SJed Brown   } else {
1148fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
1158fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
1168fe8eb27SJed Brown     *xx  = shell->right_work;
1178fe8eb27SJed Brown   }
1188fe8eb27SJed Brown   PetscFunctionReturn(0);
1198fe8eb27SJed Brown }
1208fe8eb27SJed Brown 
1218fe8eb27SJed Brown #undef __FUNCT__
1228fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleLeft"
1238fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
1248fe8eb27SJed Brown {
1258fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1268fe8eb27SJed Brown   PetscErrorCode ierr;
1278fe8eb27SJed Brown 
1288fe8eb27SJed Brown   PetscFunctionBegin;
1298fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
1308fe8eb27SJed Brown   PetscFunctionReturn(0);
1318fe8eb27SJed Brown }
1328fe8eb27SJed Brown 
1338fe8eb27SJed Brown #undef __FUNCT__
1348fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleRight"
1358fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
1368fe8eb27SJed Brown {
1378fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1388fe8eb27SJed Brown   PetscErrorCode ierr;
1398fe8eb27SJed Brown 
1408fe8eb27SJed Brown   PetscFunctionBegin;
1418fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
1428fe8eb27SJed Brown   PetscFunctionReturn(0);
1438fe8eb27SJed Brown }
1448fe8eb27SJed Brown 
1458fe8eb27SJed Brown #undef __FUNCT__
1468fe8eb27SJed Brown #define __FUNCT__ "MatShellShiftAndScale"
1478fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
1488fe8eb27SJed Brown {
1498fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1508fe8eb27SJed Brown   PetscErrorCode ierr;
1518fe8eb27SJed Brown 
1528fe8eb27SJed Brown   PetscFunctionBegin;
1538fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
1548fe8eb27SJed Brown     PetscInt          i,m;
1558fe8eb27SJed Brown     const PetscScalar *x,*d;
1568fe8eb27SJed Brown     PetscScalar       *y;
1578fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
1588fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1598fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1608fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1618fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
1628fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
1638fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1648fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
165880538faSHong Zhang   } else if (PetscAbsScalar(shell->vshift) != 0) {
1668fe8eb27SJed Brown     ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr);
167027c5fb5SJed Brown   } else if (shell->vscale != 1.0) {
168027c5fb5SJed Brown     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
1698fe8eb27SJed Brown   }
1708fe8eb27SJed Brown   PetscFunctionReturn(0);
1718fe8eb27SJed Brown }
172e51e0e81SBarry Smith 
173711e205bSSatish Balay #undef __FUNCT__
174711e205bSSatish Balay #define __FUNCT__ "MatShellGetContext"
1759d225801SJed Brown /*@
176a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
177b4fd4287SBarry Smith 
178c7fcc2eaSBarry Smith     Not Collective
179c7fcc2eaSBarry Smith 
180b4fd4287SBarry Smith     Input Parameter:
181b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
182b4fd4287SBarry Smith 
183b4fd4287SBarry Smith     Output Parameter:
184b4fd4287SBarry Smith .   ctx - the user provided context
185b4fd4287SBarry Smith 
18615091d37SBarry Smith     Level: advanced
18715091d37SBarry Smith 
188daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
189daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
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;
205ce94432eSBarry Smith   else SETERRQ(PetscObjectComm((PetscObject)mat),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);
2255edf6516SJed Brown   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
2265edf6516SJed 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;
238*e3f487b0SBarry Smith   PetscObjectState instate,outstate;
239ef66eb69SBarry Smith 
240ef66eb69SBarry Smith   PetscFunctionBegin;
2418fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
242*e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
2438fe8eb27SJed Brown   ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr);
244*e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
245*e3f487b0SBarry Smith   if (instate == outstate) {
246*e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
247*e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
248*e3f487b0SBarry Smith   }
2498fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2508fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
251ef66eb69SBarry Smith   PetscFunctionReturn(0);
252ef66eb69SBarry Smith }
253ef66eb69SBarry Smith 
254ef66eb69SBarry Smith #undef __FUNCT__
2555edf6516SJed Brown #define __FUNCT__ "MatMultAdd_Shell"
2565edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
2575edf6516SJed Brown {
2585edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
2595edf6516SJed Brown   PetscErrorCode ierr;
2605edf6516SJed Brown 
2615edf6516SJed Brown   PetscFunctionBegin;
2625edf6516SJed Brown   if (y == z) {
2635edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
2645edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
265b8430d49SBarry Smith     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
2665edf6516SJed Brown   } else {
2675edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
2685edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
2695edf6516SJed Brown   }
2705edf6516SJed Brown   PetscFunctionReturn(0);
2715edf6516SJed Brown }
2725edf6516SJed Brown 
2735edf6516SJed Brown #undef __FUNCT__
27418be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell"
27518be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
27618be62a5SSatish Balay {
27718be62a5SSatish Balay   Mat_Shell        *shell = (Mat_Shell*)A->data;
27818be62a5SSatish Balay   PetscErrorCode   ierr;
27925578ef6SJed Brown   Vec              xx;
280*e3f487b0SBarry Smith   PetscObjectState instate,outstate;
28118be62a5SSatish Balay 
28218be62a5SSatish Balay   PetscFunctionBegin;
2838fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
284*e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
2858fe8eb27SJed Brown   ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr);
286*e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
287*e3f487b0SBarry Smith   if (instate == outstate) {
288*e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
289*e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
290*e3f487b0SBarry Smith   }
2918fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
2928fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
29318be62a5SSatish Balay   PetscFunctionReturn(0);
29418be62a5SSatish Balay }
29518be62a5SSatish Balay 
29618be62a5SSatish Balay #undef __FUNCT__
2975edf6516SJed Brown #define __FUNCT__ "MatMultTransposeAdd_Shell"
2985edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
2995edf6516SJed Brown {
3005edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3015edf6516SJed Brown   PetscErrorCode ierr;
3025edf6516SJed Brown 
3035edf6516SJed Brown   PetscFunctionBegin;
3045edf6516SJed Brown   if (y == z) {
3055edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
3065edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
3075edf6516SJed Brown     ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr);
3085edf6516SJed Brown   } else {
3095edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
3105edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3115edf6516SJed Brown   }
3125edf6516SJed Brown   PetscFunctionReturn(0);
3135edf6516SJed Brown }
3145edf6516SJed Brown 
3155edf6516SJed Brown #undef __FUNCT__
31618be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell"
31718be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
31818be62a5SSatish Balay {
31918be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
32018be62a5SSatish Balay   PetscErrorCode ierr;
32118be62a5SSatish Balay 
32218be62a5SSatish Balay   PetscFunctionBegin;
32318be62a5SSatish Balay   ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr);
32418be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
3258fe8eb27SJed Brown   if (shell->dshift) {
3268fe8eb27SJed Brown     ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr);
3278fe8eb27SJed Brown   } else {
32818be62a5SSatish Balay     ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
32918be62a5SSatish Balay   }
3308fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
3318fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
33218be62a5SSatish Balay   PetscFunctionReturn(0);
33318be62a5SSatish Balay }
33418be62a5SSatish Balay 
33518be62a5SSatish Balay #undef __FUNCT__
336ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell"
337f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
338ef66eb69SBarry Smith {
339ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3408fe8eb27SJed Brown   PetscErrorCode ierr;
341b24ad042SBarry Smith 
342ef66eb69SBarry Smith   PetscFunctionBegin;
3438fe8eb27SJed Brown   if (shell->left || shell->right || shell->dshift) {
3448fe8eb27SJed Brown     if (!shell->dshift) {
3458fe8eb27SJed Brown       if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);}
3468fe8eb27SJed Brown       shell->dshift = shell->dshift_owned;
3478fe8eb27SJed Brown       ierr          = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr);
3488fe8eb27SJed Brown     } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
3498fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
3508fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
3518fe8eb27SJed Brown   } else shell->vshift += a;
3528fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
353ef66eb69SBarry Smith   PetscFunctionReturn(0);
354ef66eb69SBarry Smith }
355ef66eb69SBarry Smith 
356ef66eb69SBarry Smith #undef __FUNCT__
357ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell"
358f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
359ef66eb69SBarry Smith {
360ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3618fe8eb27SJed Brown   PetscErrorCode ierr;
362b24ad042SBarry Smith 
363ef66eb69SBarry Smith   PetscFunctionBegin;
364f4df32b1SMatthew Knepley   shell->vscale *= a;
3652205254eSKarl Rupp   if (shell->dshift) {
3662205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
3672205254eSKarl Rupp   } else shell->vshift *= a;
3688fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
3698fe8eb27SJed Brown   PetscFunctionReturn(0);
37018be62a5SSatish Balay }
3718fe8eb27SJed Brown 
3728fe8eb27SJed Brown #undef __FUNCT__
3738fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell"
3748fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
3758fe8eb27SJed Brown {
3768fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3778fe8eb27SJed Brown   PetscErrorCode ierr;
3788fe8eb27SJed Brown 
3798fe8eb27SJed Brown   PetscFunctionBegin;
3808fe8eb27SJed Brown   if (left) {
3818fe8eb27SJed Brown     if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);}
3822205254eSKarl Rupp     if (shell->left) {
3832205254eSKarl Rupp       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
3842205254eSKarl Rupp     } else {
3858fe8eb27SJed Brown       shell->left = shell->left_owned;
3868fe8eb27SJed Brown       ierr        = VecCopy(left,shell->left);CHKERRQ(ierr);
38718be62a5SSatish Balay     }
388ef66eb69SBarry Smith   }
3898fe8eb27SJed Brown   if (right) {
3908fe8eb27SJed Brown     if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);}
3912205254eSKarl Rupp     if (shell->right) {
3922205254eSKarl Rupp       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
3932205254eSKarl Rupp     } else {
3948fe8eb27SJed Brown       shell->right = shell->right_owned;
3958fe8eb27SJed Brown       ierr         = VecCopy(right,shell->right);CHKERRQ(ierr);
3968fe8eb27SJed Brown     }
3978fe8eb27SJed Brown   }
3988fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
399ef66eb69SBarry Smith   PetscFunctionReturn(0);
400ef66eb69SBarry Smith }
401ef66eb69SBarry Smith 
402ef66eb69SBarry Smith #undef __FUNCT__
403ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell"
404dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
405ef66eb69SBarry Smith {
406ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
407ef66eb69SBarry Smith 
408ef66eb69SBarry Smith   PetscFunctionBegin;
4098fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
410ef66eb69SBarry Smith     shell->vshift = 0.0;
411ef66eb69SBarry Smith     shell->vscale = 1.0;
4120298fd71SBarry Smith     shell->dshift = NULL;
4130298fd71SBarry Smith     shell->left   = NULL;
4140298fd71SBarry Smith     shell->right  = NULL;
4157fb7d96aSJed Brown     if (shell->mult) {
416ef66eb69SBarry Smith       Y->ops->mult = shell->mult;
4170298fd71SBarry Smith       shell->mult  = NULL;
4187fb7d96aSJed Brown     }
4197fb7d96aSJed Brown     if (shell->multtranspose) {
42018be62a5SSatish Balay       Y->ops->multtranspose = shell->multtranspose;
4210298fd71SBarry Smith       shell->multtranspose  = NULL;
4227fb7d96aSJed Brown     }
4237fb7d96aSJed Brown     if (shell->getdiagonal) {
42418be62a5SSatish Balay       Y->ops->getdiagonal = shell->getdiagonal;
4250298fd71SBarry Smith       shell->getdiagonal  = NULL;
4267fb7d96aSJed Brown     }
4277fb7d96aSJed Brown     shell->usingscaled = PETSC_FALSE;
428ef66eb69SBarry Smith   }
429ef66eb69SBarry Smith   PetscFunctionReturn(0);
430ef66eb69SBarry Smith }
431ef66eb69SBarry Smith 
43219fd82e9SBarry Smith extern PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
433b951964fSBarry Smith 
43409dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
43520563c6bSBarry Smith                                        0,
43620563c6bSBarry Smith                                        0,
43720563c6bSBarry Smith                                        0,
43897304618SKris Buschelman                                 /* 4*/ 0,
43920563c6bSBarry Smith                                        0,
440b951964fSBarry Smith                                        0,
441b951964fSBarry Smith                                        0,
442b951964fSBarry Smith                                        0,
443b951964fSBarry Smith                                        0,
44497304618SKris Buschelman                                 /*10*/ 0,
445b951964fSBarry Smith                                        0,
446b951964fSBarry Smith                                        0,
447b951964fSBarry Smith                                        0,
448b951964fSBarry Smith                                        0,
44997304618SKris Buschelman                                 /*15*/ 0,
450b951964fSBarry Smith                                        0,
451b951964fSBarry Smith                                        0,
4528fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
453b951964fSBarry Smith                                        0,
45497304618SKris Buschelman                                 /*20*/ 0,
455ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
456b951964fSBarry Smith                                        0,
457b951964fSBarry Smith                                        0,
458d519adbfSMatthew Knepley                                 /*24*/ 0,
459b951964fSBarry Smith                                        0,
460b951964fSBarry Smith                                        0,
461b951964fSBarry Smith                                        0,
462b951964fSBarry Smith                                        0,
463d519adbfSMatthew Knepley                                 /*29*/ 0,
464b951964fSBarry Smith                                        0,
465273d9f13SBarry Smith                                        0,
466b951964fSBarry Smith                                        0,
467b951964fSBarry Smith                                        0,
468d519adbfSMatthew Knepley                                 /*34*/ 0,
469b951964fSBarry Smith                                        0,
470b951964fSBarry Smith                                        0,
47109dc0095SBarry Smith                                        0,
47209dc0095SBarry Smith                                        0,
473d519adbfSMatthew Knepley                                 /*39*/ 0,
47409dc0095SBarry Smith                                        0,
47509dc0095SBarry Smith                                        0,
47609dc0095SBarry Smith                                        0,
47709dc0095SBarry Smith                                        0,
478d519adbfSMatthew Knepley                                 /*44*/ 0,
479ef66eb69SBarry Smith                                        MatScale_Shell,
480ef66eb69SBarry Smith                                        MatShift_Shell,
48109dc0095SBarry Smith                                        0,
48209dc0095SBarry Smith                                        0,
483f73d5cc4SBarry Smith                                 /*49*/ 0,
48409dc0095SBarry Smith                                        0,
48509dc0095SBarry Smith                                        0,
48609dc0095SBarry Smith                                        0,
48709dc0095SBarry Smith                                        0,
488d519adbfSMatthew Knepley                                 /*54*/ 0,
48909dc0095SBarry Smith                                        0,
49009dc0095SBarry Smith                                        0,
49109dc0095SBarry Smith                                        0,
49209dc0095SBarry Smith                                        0,
493d519adbfSMatthew Knepley                                 /*59*/ 0,
494b9b97703SBarry Smith                                        MatDestroy_Shell,
49509dc0095SBarry Smith                                        0,
496357abbc8SBarry Smith                                        0,
497273d9f13SBarry Smith                                        0,
498d519adbfSMatthew Knepley                                 /*64*/ 0,
499273d9f13SBarry Smith                                        0,
500273d9f13SBarry Smith                                        0,
501273d9f13SBarry Smith                                        0,
502273d9f13SBarry Smith                                        0,
503d519adbfSMatthew Knepley                                 /*69*/ 0,
504273d9f13SBarry Smith                                        0,
505c87e5d42SMatthew Knepley                                        MatConvert_Shell,
506273d9f13SBarry Smith                                        0,
50797304618SKris Buschelman                                        0,
508d519adbfSMatthew Knepley                                 /*74*/ 0,
50997304618SKris Buschelman                                        0,
51097304618SKris Buschelman                                        0,
51197304618SKris Buschelman                                        0,
51297304618SKris Buschelman                                        0,
513d519adbfSMatthew Knepley                                 /*79*/ 0,
51497304618SKris Buschelman                                        0,
51597304618SKris Buschelman                                        0,
51697304618SKris Buschelman                                        0,
517865e5f61SKris Buschelman                                        0,
518d519adbfSMatthew Knepley                                 /*84*/ 0,
519865e5f61SKris Buschelman                                        0,
520865e5f61SKris Buschelman                                        0,
521865e5f61SKris Buschelman                                        0,
522865e5f61SKris Buschelman                                        0,
523d519adbfSMatthew Knepley                                 /*89*/ 0,
524865e5f61SKris Buschelman                                        0,
525865e5f61SKris Buschelman                                        0,
526865e5f61SKris Buschelman                                        0,
527865e5f61SKris Buschelman                                        0,
528d519adbfSMatthew Knepley                                 /*94*/ 0,
529865e5f61SKris Buschelman                                        0,
530865e5f61SKris Buschelman                                        0,
5313964eb88SJed Brown                                        0,
5323964eb88SJed Brown                                        0,
5333964eb88SJed Brown                                 /*99*/ 0,
5343964eb88SJed Brown                                        0,
5353964eb88SJed Brown                                        0,
5363964eb88SJed Brown                                        0,
5373964eb88SJed Brown                                        0,
5383964eb88SJed Brown                                /*104*/ 0,
5393964eb88SJed Brown                                        0,
5403964eb88SJed Brown                                        0,
5413964eb88SJed Brown                                        0,
5423964eb88SJed Brown                                        0,
5433964eb88SJed Brown                                /*109*/ 0,
5443964eb88SJed Brown                                        0,
5453964eb88SJed Brown                                        0,
5463964eb88SJed Brown                                        0,
5473964eb88SJed Brown                                        0,
5483964eb88SJed Brown                                /*114*/ 0,
5493964eb88SJed Brown                                        0,
5503964eb88SJed Brown                                        0,
5513964eb88SJed Brown                                        0,
5523964eb88SJed Brown                                        0,
5533964eb88SJed Brown                                /*119*/ 0,
5543964eb88SJed Brown                                        0,
5553964eb88SJed Brown                                        0,
5563964eb88SJed Brown                                        0,
5573964eb88SJed Brown                                        0,
5583964eb88SJed Brown                                /*124*/ 0,
5593964eb88SJed Brown                                        0,
5603964eb88SJed Brown                                        0,
5613964eb88SJed Brown                                        0,
5623964eb88SJed Brown                                        0,
5633964eb88SJed Brown                                /*129*/ 0,
5643964eb88SJed Brown                                        0,
5653964eb88SJed Brown                                        0,
5663964eb88SJed Brown                                        0,
5673964eb88SJed Brown                                        0,
5683964eb88SJed Brown                                /*134*/ 0,
5693964eb88SJed Brown                                        0,
5703964eb88SJed Brown                                        0,
5713964eb88SJed Brown                                        0,
5723964eb88SJed Brown                                        0,
5733964eb88SJed Brown                                /*139*/ 0,
574f9426fe0SMark Adams                                        0,
5753964eb88SJed Brown                                        0
5763964eb88SJed Brown };
577273d9f13SBarry Smith 
5780bad9183SKris Buschelman /*MC
579fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
5800bad9183SKris Buschelman 
5810bad9183SKris Buschelman   Level: advanced
5820bad9183SKris Buschelman 
5830bad9183SKris Buschelman .seealso: MatCreateShell
5840bad9183SKris Buschelman M*/
5850bad9183SKris Buschelman 
586711e205bSSatish Balay #undef __FUNCT__
587711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell"
5888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
589273d9f13SBarry Smith {
590273d9f13SBarry Smith   Mat_Shell      *b;
591dfbe8321SBarry Smith   PetscErrorCode ierr;
592273d9f13SBarry Smith 
593273d9f13SBarry Smith   PetscFunctionBegin;
594273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
595273d9f13SBarry Smith 
596b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
597273d9f13SBarry Smith   A->data = (void*)b;
598273d9f13SBarry Smith 
59926283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
60026283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
601273d9f13SBarry Smith 
602273d9f13SBarry Smith   b->ctx           = 0;
603ef66eb69SBarry Smith   b->vshift        = 0.0;
604ef66eb69SBarry Smith   b->vscale        = 1.0;
605ef66eb69SBarry Smith   b->mult          = 0;
60618be62a5SSatish Balay   b->multtranspose = 0;
60718be62a5SSatish Balay   b->getdiagonal   = 0;
608273d9f13SBarry Smith   A->assembled     = PETSC_TRUE;
609210f0121SBarry Smith   A->preallocated  = PETSC_FALSE;
6102205254eSKarl Rupp 
61117667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
612273d9f13SBarry Smith   PetscFunctionReturn(0);
613273d9f13SBarry Smith }
614e51e0e81SBarry Smith 
615711e205bSSatish Balay #undef __FUNCT__
616711e205bSSatish Balay #define __FUNCT__ "MatCreateShell"
6174b828684SBarry Smith /*@C
618052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
619ff756334SLois Curfman McInnes    private data storage format.
620e51e0e81SBarry Smith 
621c7fcc2eaSBarry Smith   Collective on MPI_Comm
622c7fcc2eaSBarry Smith 
623e51e0e81SBarry Smith    Input Parameters:
624c7fcc2eaSBarry Smith +  comm - MPI communicator
625c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
626c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
627c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
628c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
629c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
630e51e0e81SBarry Smith 
631ff756334SLois Curfman McInnes    Output Parameter:
63244cd7ae7SLois Curfman McInnes .  A - the matrix
633e51e0e81SBarry Smith 
634ff2fd236SBarry Smith    Level: advanced
635ff2fd236SBarry Smith 
636f39d1f56SLois Curfman McInnes   Usage:
6377b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
638f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
639c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
640f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
641f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
642f39d1f56SLois Curfman McInnes 
643ff756334SLois Curfman McInnes    Notes:
644ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
645ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
646ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
647e51e0e81SBarry Smith 
648daf670e6SBarry Smith    Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this
649daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
650daf670e6SBarry Smith     in as the ctx argument.
651f38a8d56SBarry Smith 
652f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
653f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
654645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
655645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
656645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
657645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
658645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
659645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
660645985a0SLois Curfman McInnes    For example,
661f39d1f56SLois Curfman McInnes 
662f39d1f56SLois Curfman McInnes $
663f39d1f56SLois Curfman McInnes $     Vec x, y
6647b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
665645985a0SLois Curfman McInnes $     Mat A
666f39d1f56SLois Curfman McInnes $
667c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
668c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
669f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
670c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
671c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
672c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
673645985a0SLois Curfman McInnes $     MatMult(A,x,y);
674645985a0SLois Curfman McInnes $     MatDestroy(A);
675f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
676645985a0SLois Curfman McInnes $
677e51e0e81SBarry Smith 
6780b627109SLois Curfman McInnes .keywords: matrix, shell, create
6790b627109SLois Curfman McInnes 
680ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
681e51e0e81SBarry Smith @*/
6827087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
683e51e0e81SBarry Smith {
684dfbe8321SBarry Smith   PetscErrorCode ierr;
685ed3cc1f0SBarry Smith 
6863a40ed3dSBarry Smith   PetscFunctionBegin;
687f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
688f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
689273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
690273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
6910fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
692273d9f13SBarry Smith   PetscFunctionReturn(0);
693c7fcc2eaSBarry Smith }
694c7fcc2eaSBarry Smith 
695711e205bSSatish Balay #undef __FUNCT__
696711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext"
697c6866cfdSSatish Balay /*@
698273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
699c7fcc2eaSBarry Smith 
7003f9fe445SBarry Smith    Logically Collective on Mat
701c7fcc2eaSBarry Smith 
702273d9f13SBarry Smith     Input Parameters:
703273d9f13SBarry Smith +   mat - the shell matrix
704273d9f13SBarry Smith -   ctx - the context
705273d9f13SBarry Smith 
706273d9f13SBarry Smith    Level: advanced
707273d9f13SBarry Smith 
708daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
709daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
710273d9f13SBarry Smith 
711273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
7120bc0a719SBarry Smith @*/
7137087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
714273d9f13SBarry Smith {
715273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
716dfbe8321SBarry Smith   PetscErrorCode ierr;
717ace3abfcSBarry Smith   PetscBool      flg;
718273d9f13SBarry Smith 
719273d9f13SBarry Smith   PetscFunctionBegin;
7200700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
721251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
722273d9f13SBarry Smith   if (flg) {
723273d9f13SBarry Smith     shell->ctx = ctx;
724ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
7253a40ed3dSBarry Smith   PetscFunctionReturn(0);
726e51e0e81SBarry Smith }
727e51e0e81SBarry Smith 
728711e205bSSatish Balay #undef __FUNCT__
729711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation"
730c16cb8f2SBarry Smith /*@C
7313a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
7323a3eedf2SBarry Smith                            a shell matrix.
733e51e0e81SBarry Smith 
7343f9fe445SBarry Smith    Logically Collective on Mat
735fee21e36SBarry Smith 
736c7fcc2eaSBarry Smith     Input Parameters:
737c7fcc2eaSBarry Smith +   mat - the shell matrix
738c7fcc2eaSBarry Smith .   op - the name of the operation
739c7fcc2eaSBarry Smith -   f - the function that provides the operation.
740c7fcc2eaSBarry Smith 
74115091d37SBarry Smith    Level: advanced
74215091d37SBarry Smith 
743fae171e0SBarry Smith     Usage:
744e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
745f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
746c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
7470b627109SLois Curfman McInnes 
748a62d957aSLois Curfman McInnes     Notes:
749e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
7501c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
751a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
7521c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
753a62d957aSLois Curfman McInnes 
75425296bd5SBarry Smith     All user-provided functions (execept for MATOP_DESTROY) should have the same calling
755deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
756deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
757deebb3c3SLois Curfman McInnes     routines, e.g.,
758a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
759a62d957aSLois Curfman McInnes 
7604aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
7614aa34b0aSBarry Smith     nonzero on failure.
7624aa34b0aSBarry Smith 
763a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
764a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
765a62d957aSLois Curfman McInnes     set by MatCreateShell().
766a62d957aSLois Curfman McInnes 
7672a7a6963SBarry Smith     Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
768501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
769501d9185SBarry Smith 
770a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
771a62d957aSLois Curfman McInnes 
772ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
773e51e0e81SBarry Smith @*/
7747087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
775e51e0e81SBarry Smith {
776dfbe8321SBarry Smith   PetscErrorCode ierr;
777ace3abfcSBarry Smith   PetscBool      flg;
778273d9f13SBarry Smith 
7793a40ed3dSBarry Smith   PetscFunctionBegin;
7800700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
7815edf6516SJed Brown   switch (op) {
7825edf6516SJed Brown   case MATOP_DESTROY:
783251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
784273d9f13SBarry Smith     if (flg) {
785a62d957aSLois Curfman McInnes       Mat_Shell *shell = (Mat_Shell*)mat->data;
7866849ba73SBarry Smith       shell->destroy = (PetscErrorCode (*)(Mat))f;
7876849ba73SBarry Smith     } else mat->ops->destroy = (PetscErrorCode (*)(Mat))f;
7885edf6516SJed Brown     break;
7895edf6516SJed Brown   case MATOP_VIEW:
7905edf6516SJed Brown     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
7915edf6516SJed Brown     break;
7925edf6516SJed Brown   case MATOP_MULT:
7935edf6516SJed Brown     mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
7945edf6516SJed Brown     if (!mat->ops->multadd) mat->ops->multadd = MatMultAdd_Shell;
7955edf6516SJed Brown     break;
7965edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
7975edf6516SJed Brown     mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
7985edf6516SJed Brown     if (!mat->ops->multtransposeadd) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell;
7995edf6516SJed Brown     break;
8005edf6516SJed Brown   default:
8015edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
802a62d957aSLois Curfman McInnes   }
8033a40ed3dSBarry Smith   PetscFunctionReturn(0);
804e51e0e81SBarry Smith }
805f0479e8cSBarry Smith 
806711e205bSSatish Balay #undef __FUNCT__
807711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation"
808d4bb536fSBarry Smith /*@C
809d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
810d4bb536fSBarry Smith 
811c7fcc2eaSBarry Smith     Not Collective
812c7fcc2eaSBarry Smith 
813d4bb536fSBarry Smith     Input Parameters:
814c7fcc2eaSBarry Smith +   mat - the shell matrix
815c7fcc2eaSBarry Smith -   op - the name of the operation
816d4bb536fSBarry Smith 
817d4bb536fSBarry Smith     Output Parameter:
818d4bb536fSBarry Smith .   f - the function that provides the operation.
819d4bb536fSBarry Smith 
82015091d37SBarry Smith     Level: advanced
82115091d37SBarry Smith 
822d4bb536fSBarry Smith     Notes:
823e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
824d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
825d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
826d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
827d4bb536fSBarry Smith 
828d4bb536fSBarry Smith     All user-provided functions have the same calling
829d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
830d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
831d4bb536fSBarry Smith     routines, e.g.,
832d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
833d4bb536fSBarry Smith 
834d4bb536fSBarry Smith     Within each user-defined routine, the user should call
835d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
836d4bb536fSBarry Smith     set by MatCreateShell().
837d4bb536fSBarry Smith 
838d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
839d4bb536fSBarry Smith 
840ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
841d4bb536fSBarry Smith @*/
8427087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
843d4bb536fSBarry Smith {
844dfbe8321SBarry Smith   PetscErrorCode ierr;
845ace3abfcSBarry Smith   PetscBool      flg;
846273d9f13SBarry Smith 
8473a40ed3dSBarry Smith   PetscFunctionBegin;
8480700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
849d4bb536fSBarry Smith   if (op == MATOP_DESTROY) {
850251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
851273d9f13SBarry Smith     if (flg) {
852d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
853c134de8dSSatish Balay       *f = (void (*)(void))shell->destroy;
854c7fcc2eaSBarry Smith     } else {
855c134de8dSSatish Balay       *f = (void (*)(void))mat->ops->destroy;
856d4bb536fSBarry Smith     }
857c7fcc2eaSBarry Smith   } else if (op == MATOP_VIEW) {
858c134de8dSSatish Balay     *f = (void (*)(void))mat->ops->view;
859c7fcc2eaSBarry Smith   } else {
860c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
861d4bb536fSBarry Smith   }
8623a40ed3dSBarry Smith   PetscFunctionReturn(0);
863d4bb536fSBarry Smith }
864d4bb536fSBarry Smith 
865