xref: /petsc/src/mat/impls/shell/shell.c (revision cc2e6a90c05b27ffec69cb207fe793d447f14420)
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;
238e3f487b0SBarry Smith   PetscObjectState instate,outstate;
239ef66eb69SBarry Smith 
240ef66eb69SBarry Smith   PetscFunctionBegin;
2418fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
242e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
2438fe8eb27SJed Brown   ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr);
244e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
245e3f487b0SBarry Smith   if (instate == outstate) {
246e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
247e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
248e3f487b0SBarry 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;
280e3f487b0SBarry Smith   PetscObjectState instate,outstate;
28118be62a5SSatish Balay 
28218be62a5SSatish Balay   PetscFunctionBegin;
2838fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
284e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
2858fe8eb27SJed Brown   ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr);
286e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
287e3f487b0SBarry Smith   if (instate == outstate) {
288e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
289e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
290e3f487b0SBarry 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 
432*cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
433b951964fSBarry Smith 
4343b49f96aSBarry Smith #undef __FUNCT__
4353b49f96aSBarry Smith #define __FUNCT__ "MatMissingDiagonal_Shell"
4363b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
4373b49f96aSBarry Smith {
4383b49f96aSBarry Smith   PetscFunctionBegin;
4393b49f96aSBarry Smith   *missing = PETSC_FALSE;
4403b49f96aSBarry Smith   PetscFunctionReturn(0);
4413b49f96aSBarry Smith }
4423b49f96aSBarry Smith 
44309dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
44420563c6bSBarry Smith                                        0,
44520563c6bSBarry Smith                                        0,
44620563c6bSBarry Smith                                        0,
44797304618SKris Buschelman                                 /* 4*/ 0,
44820563c6bSBarry Smith                                        0,
449b951964fSBarry Smith                                        0,
450b951964fSBarry Smith                                        0,
451b951964fSBarry Smith                                        0,
452b951964fSBarry Smith                                        0,
45397304618SKris Buschelman                                 /*10*/ 0,
454b951964fSBarry Smith                                        0,
455b951964fSBarry Smith                                        0,
456b951964fSBarry Smith                                        0,
457b951964fSBarry Smith                                        0,
45897304618SKris Buschelman                                 /*15*/ 0,
459b951964fSBarry Smith                                        0,
460b951964fSBarry Smith                                        0,
4618fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
462b951964fSBarry Smith                                        0,
46397304618SKris Buschelman                                 /*20*/ 0,
464ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
465b951964fSBarry Smith                                        0,
466b951964fSBarry Smith                                        0,
467d519adbfSMatthew Knepley                                 /*24*/ 0,
468b951964fSBarry Smith                                        0,
469b951964fSBarry Smith                                        0,
470b951964fSBarry Smith                                        0,
471b951964fSBarry Smith                                        0,
472d519adbfSMatthew Knepley                                 /*29*/ 0,
473b951964fSBarry Smith                                        0,
474273d9f13SBarry Smith                                        0,
475b951964fSBarry Smith                                        0,
476b951964fSBarry Smith                                        0,
477d519adbfSMatthew Knepley                                 /*34*/ 0,
478b951964fSBarry Smith                                        0,
479b951964fSBarry Smith                                        0,
48009dc0095SBarry Smith                                        0,
48109dc0095SBarry Smith                                        0,
482d519adbfSMatthew Knepley                                 /*39*/ 0,
48309dc0095SBarry Smith                                        0,
48409dc0095SBarry Smith                                        0,
48509dc0095SBarry Smith                                        0,
48609dc0095SBarry Smith                                        0,
487d519adbfSMatthew Knepley                                 /*44*/ 0,
488ef66eb69SBarry Smith                                        MatScale_Shell,
489ef66eb69SBarry Smith                                        MatShift_Shell,
49009dc0095SBarry Smith                                        0,
49109dc0095SBarry Smith                                        0,
492f73d5cc4SBarry Smith                                 /*49*/ 0,
49309dc0095SBarry Smith                                        0,
49409dc0095SBarry Smith                                        0,
49509dc0095SBarry Smith                                        0,
49609dc0095SBarry Smith                                        0,
497d519adbfSMatthew Knepley                                 /*54*/ 0,
49809dc0095SBarry Smith                                        0,
49909dc0095SBarry Smith                                        0,
50009dc0095SBarry Smith                                        0,
50109dc0095SBarry Smith                                        0,
502d519adbfSMatthew Knepley                                 /*59*/ 0,
503b9b97703SBarry Smith                                        MatDestroy_Shell,
50409dc0095SBarry Smith                                        0,
505357abbc8SBarry Smith                                        0,
506273d9f13SBarry Smith                                        0,
507d519adbfSMatthew Knepley                                 /*64*/ 0,
508273d9f13SBarry Smith                                        0,
509273d9f13SBarry Smith                                        0,
510273d9f13SBarry Smith                                        0,
511273d9f13SBarry Smith                                        0,
512d519adbfSMatthew Knepley                                 /*69*/ 0,
513273d9f13SBarry Smith                                        0,
514c87e5d42SMatthew Knepley                                        MatConvert_Shell,
515273d9f13SBarry Smith                                        0,
51697304618SKris Buschelman                                        0,
517d519adbfSMatthew Knepley                                 /*74*/ 0,
51897304618SKris Buschelman                                        0,
51997304618SKris Buschelman                                        0,
52097304618SKris Buschelman                                        0,
52197304618SKris Buschelman                                        0,
522d519adbfSMatthew Knepley                                 /*79*/ 0,
52397304618SKris Buschelman                                        0,
52497304618SKris Buschelman                                        0,
52597304618SKris Buschelman                                        0,
526865e5f61SKris Buschelman                                        0,
527d519adbfSMatthew Knepley                                 /*84*/ 0,
528865e5f61SKris Buschelman                                        0,
529865e5f61SKris Buschelman                                        0,
530865e5f61SKris Buschelman                                        0,
531865e5f61SKris Buschelman                                        0,
532d519adbfSMatthew Knepley                                 /*89*/ 0,
533865e5f61SKris Buschelman                                        0,
534865e5f61SKris Buschelman                                        0,
535865e5f61SKris Buschelman                                        0,
536865e5f61SKris Buschelman                                        0,
537d519adbfSMatthew Knepley                                 /*94*/ 0,
538865e5f61SKris Buschelman                                        0,
539865e5f61SKris Buschelman                                        0,
5403964eb88SJed Brown                                        0,
5413964eb88SJed Brown                                        0,
5423964eb88SJed Brown                                 /*99*/ 0,
5433964eb88SJed Brown                                        0,
5443964eb88SJed Brown                                        0,
5453964eb88SJed Brown                                        0,
5463964eb88SJed Brown                                        0,
5473964eb88SJed Brown                                /*104*/ 0,
5483964eb88SJed Brown                                        0,
5493964eb88SJed Brown                                        0,
5503964eb88SJed Brown                                        0,
5513964eb88SJed Brown                                        0,
5523964eb88SJed Brown                                /*109*/ 0,
5533964eb88SJed Brown                                        0,
5543964eb88SJed Brown                                        0,
5553964eb88SJed Brown                                        0,
5563b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
5573964eb88SJed Brown                                /*114*/ 0,
5583964eb88SJed Brown                                        0,
5593964eb88SJed Brown                                        0,
5603964eb88SJed Brown                                        0,
5613964eb88SJed Brown                                        0,
5623964eb88SJed Brown                                /*119*/ 0,
5633964eb88SJed Brown                                        0,
5643964eb88SJed Brown                                        0,
5653964eb88SJed Brown                                        0,
5663964eb88SJed Brown                                        0,
5673964eb88SJed Brown                                /*124*/ 0,
5683964eb88SJed Brown                                        0,
5693964eb88SJed Brown                                        0,
5703964eb88SJed Brown                                        0,
5713964eb88SJed Brown                                        0,
5723964eb88SJed Brown                                /*129*/ 0,
5733964eb88SJed Brown                                        0,
5743964eb88SJed Brown                                        0,
5753964eb88SJed Brown                                        0,
5763964eb88SJed Brown                                        0,
5773964eb88SJed Brown                                /*134*/ 0,
5783964eb88SJed Brown                                        0,
5793964eb88SJed Brown                                        0,
5803964eb88SJed Brown                                        0,
5813964eb88SJed Brown                                        0,
5823964eb88SJed Brown                                /*139*/ 0,
583f9426fe0SMark Adams                                        0,
5843964eb88SJed Brown                                        0
5853964eb88SJed Brown };
586273d9f13SBarry Smith 
5870bad9183SKris Buschelman /*MC
588fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
5890bad9183SKris Buschelman 
5900bad9183SKris Buschelman   Level: advanced
5910bad9183SKris Buschelman 
5920bad9183SKris Buschelman .seealso: MatCreateShell
5930bad9183SKris Buschelman M*/
5940bad9183SKris Buschelman 
595711e205bSSatish Balay #undef __FUNCT__
596711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell"
5978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
598273d9f13SBarry Smith {
599273d9f13SBarry Smith   Mat_Shell      *b;
600dfbe8321SBarry Smith   PetscErrorCode ierr;
601273d9f13SBarry Smith 
602273d9f13SBarry Smith   PetscFunctionBegin;
603273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
604273d9f13SBarry Smith 
605b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
606273d9f13SBarry Smith   A->data = (void*)b;
607273d9f13SBarry Smith 
60826283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
60926283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
610273d9f13SBarry Smith 
611273d9f13SBarry Smith   b->ctx           = 0;
612ef66eb69SBarry Smith   b->vshift        = 0.0;
613ef66eb69SBarry Smith   b->vscale        = 1.0;
614ef66eb69SBarry Smith   b->mult          = 0;
61518be62a5SSatish Balay   b->multtranspose = 0;
61618be62a5SSatish Balay   b->getdiagonal   = 0;
617273d9f13SBarry Smith   A->assembled     = PETSC_TRUE;
618210f0121SBarry Smith   A->preallocated  = PETSC_FALSE;
6192205254eSKarl Rupp 
62017667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
621273d9f13SBarry Smith   PetscFunctionReturn(0);
622273d9f13SBarry Smith }
623e51e0e81SBarry Smith 
624711e205bSSatish Balay #undef __FUNCT__
625711e205bSSatish Balay #define __FUNCT__ "MatCreateShell"
6264b828684SBarry Smith /*@C
627052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
628ff756334SLois Curfman McInnes    private data storage format.
629e51e0e81SBarry Smith 
630c7fcc2eaSBarry Smith   Collective on MPI_Comm
631c7fcc2eaSBarry Smith 
632e51e0e81SBarry Smith    Input Parameters:
633c7fcc2eaSBarry Smith +  comm - MPI communicator
634c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
635c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
636c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
637c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
638c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
639e51e0e81SBarry Smith 
640ff756334SLois Curfman McInnes    Output Parameter:
64144cd7ae7SLois Curfman McInnes .  A - the matrix
642e51e0e81SBarry Smith 
643ff2fd236SBarry Smith    Level: advanced
644ff2fd236SBarry Smith 
645f39d1f56SLois Curfman McInnes   Usage:
6467b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
647f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
648c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
649f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
650f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
651f39d1f56SLois Curfman McInnes 
652ff756334SLois Curfman McInnes    Notes:
653ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
654ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
655ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
656e51e0e81SBarry Smith 
657daf670e6SBarry Smith    Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this
658daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
659daf670e6SBarry Smith     in as the ctx argument.
660f38a8d56SBarry Smith 
661f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
662f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
663645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
664645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
665645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
666645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
667645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
668645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
669645985a0SLois Curfman McInnes    For example,
670f39d1f56SLois Curfman McInnes 
671f39d1f56SLois Curfman McInnes $
672f39d1f56SLois Curfman McInnes $     Vec x, y
6737b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
674645985a0SLois Curfman McInnes $     Mat A
675f39d1f56SLois Curfman McInnes $
676c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
677c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
678f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
679c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
680c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
681c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
682645985a0SLois Curfman McInnes $     MatMult(A,x,y);
683645985a0SLois Curfman McInnes $     MatDestroy(A);
684f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
685645985a0SLois Curfman McInnes $
686e51e0e81SBarry Smith 
6870b627109SLois Curfman McInnes .keywords: matrix, shell, create
6880b627109SLois Curfman McInnes 
689ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
690e51e0e81SBarry Smith @*/
6917087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
692e51e0e81SBarry Smith {
693dfbe8321SBarry Smith   PetscErrorCode ierr;
694ed3cc1f0SBarry Smith 
6953a40ed3dSBarry Smith   PetscFunctionBegin;
696f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
697f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
698273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
699273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
7000fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
701273d9f13SBarry Smith   PetscFunctionReturn(0);
702c7fcc2eaSBarry Smith }
703c7fcc2eaSBarry Smith 
704711e205bSSatish Balay #undef __FUNCT__
705711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext"
706c6866cfdSSatish Balay /*@
707273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
708c7fcc2eaSBarry Smith 
7093f9fe445SBarry Smith    Logically Collective on Mat
710c7fcc2eaSBarry Smith 
711273d9f13SBarry Smith     Input Parameters:
712273d9f13SBarry Smith +   mat - the shell matrix
713273d9f13SBarry Smith -   ctx - the context
714273d9f13SBarry Smith 
715273d9f13SBarry Smith    Level: advanced
716273d9f13SBarry Smith 
717daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
718daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
719273d9f13SBarry Smith 
720273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
7210bc0a719SBarry Smith @*/
7227087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
723273d9f13SBarry Smith {
724273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
725dfbe8321SBarry Smith   PetscErrorCode ierr;
726ace3abfcSBarry Smith   PetscBool      flg;
727273d9f13SBarry Smith 
728273d9f13SBarry Smith   PetscFunctionBegin;
7290700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
730251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
731273d9f13SBarry Smith   if (flg) {
732273d9f13SBarry Smith     shell->ctx = ctx;
733ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
7343a40ed3dSBarry Smith   PetscFunctionReturn(0);
735e51e0e81SBarry Smith }
736e51e0e81SBarry Smith 
737711e205bSSatish Balay #undef __FUNCT__
738711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation"
739c16cb8f2SBarry Smith /*@C
7403a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
7413a3eedf2SBarry Smith                            a shell matrix.
742e51e0e81SBarry Smith 
7433f9fe445SBarry Smith    Logically Collective on Mat
744fee21e36SBarry Smith 
745c7fcc2eaSBarry Smith     Input Parameters:
746c7fcc2eaSBarry Smith +   mat - the shell matrix
747c7fcc2eaSBarry Smith .   op - the name of the operation
748c7fcc2eaSBarry Smith -   f - the function that provides the operation.
749c7fcc2eaSBarry Smith 
75015091d37SBarry Smith    Level: advanced
75115091d37SBarry Smith 
752fae171e0SBarry Smith     Usage:
753e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
754f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
755c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
7560b627109SLois Curfman McInnes 
757a62d957aSLois Curfman McInnes     Notes:
758e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
7591c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
760a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
7611c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
762a62d957aSLois Curfman McInnes 
76325296bd5SBarry Smith     All user-provided functions (execept for MATOP_DESTROY) should have the same calling
764deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
765deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
766deebb3c3SLois Curfman McInnes     routines, e.g.,
767a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
768a62d957aSLois Curfman McInnes 
7694aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
7704aa34b0aSBarry Smith     nonzero on failure.
7714aa34b0aSBarry Smith 
772a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
773a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
774a62d957aSLois Curfman McInnes     set by MatCreateShell().
775a62d957aSLois Curfman McInnes 
7762a7a6963SBarry Smith     Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
777501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
778501d9185SBarry Smith 
779a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
780a62d957aSLois Curfman McInnes 
781ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
782e51e0e81SBarry Smith @*/
7837087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
784e51e0e81SBarry Smith {
785dfbe8321SBarry Smith   PetscErrorCode ierr;
786ace3abfcSBarry Smith   PetscBool      flg;
787273d9f13SBarry Smith 
7883a40ed3dSBarry Smith   PetscFunctionBegin;
7890700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
7905edf6516SJed Brown   switch (op) {
7915edf6516SJed Brown   case MATOP_DESTROY:
792251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
793273d9f13SBarry Smith     if (flg) {
794a62d957aSLois Curfman McInnes       Mat_Shell *shell = (Mat_Shell*)mat->data;
7956849ba73SBarry Smith       shell->destroy = (PetscErrorCode (*)(Mat))f;
7966849ba73SBarry Smith     } else mat->ops->destroy = (PetscErrorCode (*)(Mat))f;
7975edf6516SJed Brown     break;
7985edf6516SJed Brown   case MATOP_VIEW:
7995edf6516SJed Brown     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
8005edf6516SJed Brown     break;
8015edf6516SJed Brown   case MATOP_MULT:
8025edf6516SJed Brown     mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
8035edf6516SJed Brown     if (!mat->ops->multadd) mat->ops->multadd = MatMultAdd_Shell;
8045edf6516SJed Brown     break;
8055edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
8065edf6516SJed Brown     mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
8075edf6516SJed Brown     if (!mat->ops->multtransposeadd) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell;
8085edf6516SJed Brown     break;
8095edf6516SJed Brown   default:
8105edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
811a62d957aSLois Curfman McInnes   }
8123a40ed3dSBarry Smith   PetscFunctionReturn(0);
813e51e0e81SBarry Smith }
814f0479e8cSBarry Smith 
815711e205bSSatish Balay #undef __FUNCT__
816711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation"
817d4bb536fSBarry Smith /*@C
818d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
819d4bb536fSBarry Smith 
820c7fcc2eaSBarry Smith     Not Collective
821c7fcc2eaSBarry Smith 
822d4bb536fSBarry Smith     Input Parameters:
823c7fcc2eaSBarry Smith +   mat - the shell matrix
824c7fcc2eaSBarry Smith -   op - the name of the operation
825d4bb536fSBarry Smith 
826d4bb536fSBarry Smith     Output Parameter:
827d4bb536fSBarry Smith .   f - the function that provides the operation.
828d4bb536fSBarry Smith 
82915091d37SBarry Smith     Level: advanced
83015091d37SBarry Smith 
831d4bb536fSBarry Smith     Notes:
832e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
833d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
834d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
835d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
836d4bb536fSBarry Smith 
837d4bb536fSBarry Smith     All user-provided functions have the same calling
838d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
839d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
840d4bb536fSBarry Smith     routines, e.g.,
841d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
842d4bb536fSBarry Smith 
843d4bb536fSBarry Smith     Within each user-defined routine, the user should call
844d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
845d4bb536fSBarry Smith     set by MatCreateShell().
846d4bb536fSBarry Smith 
847d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
848d4bb536fSBarry Smith 
849ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
850d4bb536fSBarry Smith @*/
8517087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
852d4bb536fSBarry Smith {
853dfbe8321SBarry Smith   PetscErrorCode ierr;
854ace3abfcSBarry Smith   PetscBool      flg;
855273d9f13SBarry Smith 
8563a40ed3dSBarry Smith   PetscFunctionBegin;
8570700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
858d4bb536fSBarry Smith   if (op == MATOP_DESTROY) {
859251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
860273d9f13SBarry Smith     if (flg) {
861d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
862c134de8dSSatish Balay       *f = (void (*)(void))shell->destroy;
863c7fcc2eaSBarry Smith     } else {
864c134de8dSSatish Balay       *f = (void (*)(void))mat->ops->destroy;
865d4bb536fSBarry Smith     }
866c7fcc2eaSBarry Smith   } else if (op == MATOP_VIEW) {
867c134de8dSSatish Balay     *f = (void (*)(void))mat->ops->view;
868c7fcc2eaSBarry Smith   } else {
869c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
870d4bb536fSBarry Smith   }
8713a40ed3dSBarry Smith   PetscFunctionReturn(0);
872d4bb536fSBarry Smith }
873d4bb536fSBarry Smith 
874