xref: /petsc/src/mat/impls/shell/shell.c (revision 7fabda3f22b6e7b364ea56fd6b7f4ad1c93c0a8b)
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*/
9e51e0e81SBarry Smith 
1020563c6bSBarry Smith typedef struct {
116849ba73SBarry Smith   PetscErrorCode (*destroy)(Mat);
126849ba73SBarry Smith   PetscErrorCode (*mult)(Mat,Vec,Vec);
1318be62a5SSatish Balay   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
1418be62a5SSatish Balay   PetscErrorCode (*getdiagonal)(Mat,Vec);
156464bf51SAlex Fikl   PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode);
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__
232*7fabda3fSAlex Fikl #define __FUNCT__ "MatCopy_Shell"
233*7fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
234*7fabda3fSAlex Fikl {
235*7fabda3fSAlex Fikl   Mat_Shell       *shellA,*shellB;
236*7fabda3fSAlex Fikl   PetscBool       flg;
237*7fabda3fSAlex Fikl   PetscErrorCode  ierr;
238*7fabda3fSAlex Fikl 
239*7fabda3fSAlex Fikl   PetscFunctionBegin;
240*7fabda3fSAlex Fikl   ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&flg);CHKERRQ(ierr);
241*7fabda3fSAlex Fikl   if(!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL");
242*7fabda3fSAlex Fikl 
243*7fabda3fSAlex Fikl   shellA = (Mat_Shell*)A->data;
244*7fabda3fSAlex Fikl   shellB = (Mat_Shell*)B->data;
245*7fabda3fSAlex Fikl   if (shellA->usingscaled) {
246*7fabda3fSAlex Fikl     ierr = MatShellUseScaledMethods(B);CHKERRQ(ierr);
247*7fabda3fSAlex Fikl 
248*7fabda3fSAlex Fikl     shellB->vscale = shellA->vscale;
249*7fabda3fSAlex Fikl     shellB->vshift = shellA->vshift;
250*7fabda3fSAlex Fikl     if (shellA->dshift_owned) {
251*7fabda3fSAlex Fikl       if (!shellB->dshift_owned) {
252*7fabda3fSAlex Fikl         ierr = VecDuplicate(shellA->dshift_owned,&shellB->dshift_owned);CHKERRQ(ierr);
253*7fabda3fSAlex Fikl       }
254*7fabda3fSAlex Fikl       ierr = VecCopy(shellA->dshift_owned,shellB->dshift_owned);CHKERRQ(ierr);
255*7fabda3fSAlex Fikl       shellB->dshift = shellB->dshift_owned;
256*7fabda3fSAlex Fikl     } else {
257*7fabda3fSAlex Fikl       shellB->dshift = NULL;
258*7fabda3fSAlex Fikl     }
259*7fabda3fSAlex Fikl     if (shellA->left_owned) {
260*7fabda3fSAlex Fikl       if (!shellB->left_owned) {
261*7fabda3fSAlex Fikl         ierr = VecDuplicate(shellA->left_owned,&shellB->left_owned);CHKERRQ(ierr);
262*7fabda3fSAlex Fikl       }
263*7fabda3fSAlex Fikl       ierr = VecCopy(shellA->left_owned,shellB->left_owned);CHKERRQ(ierr);
264*7fabda3fSAlex Fikl       shellB->left = shellB->left_owned;
265*7fabda3fSAlex Fikl     } else {
266*7fabda3fSAlex Fikl       shellB->left = NULL;
267*7fabda3fSAlex Fikl     }
268*7fabda3fSAlex Fikl     if (shellA->right_owned) {
269*7fabda3fSAlex Fikl       if (!shellB->right_owned) {
270*7fabda3fSAlex Fikl         ierr = VecDuplicate(shellA->right_owned,&shellB->right_owned);CHKERRQ(ierr);
271*7fabda3fSAlex Fikl       }
272*7fabda3fSAlex Fikl       ierr = VecCopy(shellA->right_owned,shellB->right_owned);CHKERRQ(ierr);
273*7fabda3fSAlex Fikl       shellB->right = shellB->right_owned;
274*7fabda3fSAlex Fikl     } else {
275*7fabda3fSAlex Fikl       shellB->right = NULL;
276*7fabda3fSAlex Fikl     }
277*7fabda3fSAlex Fikl   }
278*7fabda3fSAlex Fikl   PetscFunctionReturn(0);
279*7fabda3fSAlex Fikl }
280*7fabda3fSAlex Fikl 
281*7fabda3fSAlex Fikl #undef __FUNCT__
282ef66eb69SBarry Smith #define __FUNCT__ "MatMult_Shell"
283dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
284ef66eb69SBarry Smith {
285ef66eb69SBarry Smith   Mat_Shell        *shell = (Mat_Shell*)A->data;
286dfbe8321SBarry Smith   PetscErrorCode   ierr;
28725578ef6SJed Brown   Vec              xx;
288e3f487b0SBarry Smith   PetscObjectState instate,outstate;
289ef66eb69SBarry Smith 
290ef66eb69SBarry Smith   PetscFunctionBegin;
2918fe8eb27SJed Brown   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
292e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
2938fe8eb27SJed Brown   ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr);
294e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
295e3f487b0SBarry Smith   if (instate == outstate) {
296e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
297e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
298e3f487b0SBarry Smith   }
2998fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
3008fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
301ef66eb69SBarry Smith   PetscFunctionReturn(0);
302ef66eb69SBarry Smith }
303ef66eb69SBarry Smith 
304ef66eb69SBarry Smith #undef __FUNCT__
3055edf6516SJed Brown #define __FUNCT__ "MatMultAdd_Shell"
3065edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3075edf6516SJed Brown {
3085edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3095edf6516SJed Brown   PetscErrorCode ierr;
3105edf6516SJed Brown 
3115edf6516SJed Brown   PetscFunctionBegin;
3125edf6516SJed Brown   if (y == z) {
3135edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
3145edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
315b8430d49SBarry Smith     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
3165edf6516SJed Brown   } else {
3175edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
3185edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3195edf6516SJed Brown   }
3205edf6516SJed Brown   PetscFunctionReturn(0);
3215edf6516SJed Brown }
3225edf6516SJed Brown 
3235edf6516SJed Brown #undef __FUNCT__
32418be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell"
32518be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
32618be62a5SSatish Balay {
32718be62a5SSatish Balay   Mat_Shell        *shell = (Mat_Shell*)A->data;
32818be62a5SSatish Balay   PetscErrorCode   ierr;
32925578ef6SJed Brown   Vec              xx;
330e3f487b0SBarry Smith   PetscObjectState instate,outstate;
33118be62a5SSatish Balay 
33218be62a5SSatish Balay   PetscFunctionBegin;
3338fe8eb27SJed Brown   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
334e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
3358fe8eb27SJed Brown   ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr);
336e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
337e3f487b0SBarry Smith   if (instate == outstate) {
338e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
339e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
340e3f487b0SBarry Smith   }
3418fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
3428fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
34318be62a5SSatish Balay   PetscFunctionReturn(0);
34418be62a5SSatish Balay }
34518be62a5SSatish Balay 
34618be62a5SSatish Balay #undef __FUNCT__
3475edf6516SJed Brown #define __FUNCT__ "MatMultTransposeAdd_Shell"
3485edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
3495edf6516SJed Brown {
3505edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
3515edf6516SJed Brown   PetscErrorCode ierr;
3525edf6516SJed Brown 
3535edf6516SJed Brown   PetscFunctionBegin;
3545edf6516SJed Brown   if (y == z) {
3555edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
3565edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
3575edf6516SJed Brown     ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr);
3585edf6516SJed Brown   } else {
3595edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
3605edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
3615edf6516SJed Brown   }
3625edf6516SJed Brown   PetscFunctionReturn(0);
3635edf6516SJed Brown }
3645edf6516SJed Brown 
3655edf6516SJed Brown #undef __FUNCT__
36618be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell"
36718be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
36818be62a5SSatish Balay {
36918be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
37018be62a5SSatish Balay   PetscErrorCode ierr;
37118be62a5SSatish Balay 
37218be62a5SSatish Balay   PetscFunctionBegin;
37318be62a5SSatish Balay   ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr);
37418be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
3758fe8eb27SJed Brown   if (shell->dshift) {
3768fe8eb27SJed Brown     ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr);
3778fe8eb27SJed Brown   } else {
37818be62a5SSatish Balay     ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
37918be62a5SSatish Balay   }
3808fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
3818fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
38218be62a5SSatish Balay   PetscFunctionReturn(0);
38318be62a5SSatish Balay }
38418be62a5SSatish Balay 
38518be62a5SSatish Balay #undef __FUNCT__
386ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell"
387f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
388ef66eb69SBarry Smith {
389ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
3908fe8eb27SJed Brown   PetscErrorCode ierr;
391b24ad042SBarry Smith 
392ef66eb69SBarry Smith   PetscFunctionBegin;
3938fe8eb27SJed Brown   if (shell->left || shell->right || shell->dshift) {
3948fe8eb27SJed Brown     if (!shell->dshift) {
3958fe8eb27SJed Brown       if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);}
3968fe8eb27SJed Brown       shell->dshift = shell->dshift_owned;
3978fe8eb27SJed Brown       ierr          = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr);
3988fe8eb27SJed Brown     } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
3998fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
4008fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
4018fe8eb27SJed Brown   } else shell->vshift += a;
4028fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
403ef66eb69SBarry Smith   PetscFunctionReturn(0);
404ef66eb69SBarry Smith }
405ef66eb69SBarry Smith 
406ef66eb69SBarry Smith #undef __FUNCT__
4076464bf51SAlex Fikl #define __FUNCT__ "MatDiagonalSet_Shell"
4086464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
4096464bf51SAlex Fikl {
4106464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
4116464bf51SAlex Fikl   PetscErrorCode ierr;
4126464bf51SAlex Fikl 
4136464bf51SAlex Fikl   PetscFunctionBegin;
4146464bf51SAlex Fikl   if (shell->vscale != 1.0 || shell->left || shell->right) {
4156464bf51SAlex Fikl     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported on scaled matrix");
4166464bf51SAlex Fikl   }
4176464bf51SAlex Fikl 
4186464bf51SAlex Fikl   if (shell->diagonalset) { ierr = (*shell->diagonalset)(A,D,ins);CHKERRQ(ierr); }
4196464bf51SAlex Fikl   shell->vshift = 0.0;
4206464bf51SAlex Fikl   if (shell->dshift) { ierr = VecZeroEntries(shell->dshift);CHKERRQ(ierr); }
4216464bf51SAlex Fikl   PetscFunctionReturn(0);
4226464bf51SAlex Fikl }
4236464bf51SAlex Fikl 
4246464bf51SAlex Fikl #undef __FUNCT__
425ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell"
426f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
427ef66eb69SBarry Smith {
428ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4298fe8eb27SJed Brown   PetscErrorCode ierr;
430b24ad042SBarry Smith 
431ef66eb69SBarry Smith   PetscFunctionBegin;
432f4df32b1SMatthew Knepley   shell->vscale *= a;
4332205254eSKarl Rupp   if (shell->dshift) {
4342205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
4352205254eSKarl Rupp   } else shell->vshift *= a;
4368fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
4378fe8eb27SJed Brown   PetscFunctionReturn(0);
43818be62a5SSatish Balay }
4398fe8eb27SJed Brown 
4408fe8eb27SJed Brown #undef __FUNCT__
4418fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell"
4428fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
4438fe8eb27SJed Brown {
4448fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
4458fe8eb27SJed Brown   PetscErrorCode ierr;
4468fe8eb27SJed Brown 
4478fe8eb27SJed Brown   PetscFunctionBegin;
4488fe8eb27SJed Brown   if (left) {
4498fe8eb27SJed Brown     if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);}
4502205254eSKarl Rupp     if (shell->left) {
4512205254eSKarl Rupp       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
4522205254eSKarl Rupp     } else {
4538fe8eb27SJed Brown       shell->left = shell->left_owned;
4548fe8eb27SJed Brown       ierr        = VecCopy(left,shell->left);CHKERRQ(ierr);
45518be62a5SSatish Balay     }
456ef66eb69SBarry Smith   }
4578fe8eb27SJed Brown   if (right) {
4588fe8eb27SJed Brown     if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);}
4592205254eSKarl Rupp     if (shell->right) {
4602205254eSKarl Rupp       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
4612205254eSKarl Rupp     } else {
4628fe8eb27SJed Brown       shell->right = shell->right_owned;
4638fe8eb27SJed Brown       ierr         = VecCopy(right,shell->right);CHKERRQ(ierr);
4648fe8eb27SJed Brown     }
4658fe8eb27SJed Brown   }
4668fe8eb27SJed Brown   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
467ef66eb69SBarry Smith   PetscFunctionReturn(0);
468ef66eb69SBarry Smith }
469ef66eb69SBarry Smith 
470ef66eb69SBarry Smith #undef __FUNCT__
471ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell"
472dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
473ef66eb69SBarry Smith {
474ef66eb69SBarry Smith   Mat_Shell *shell = (Mat_Shell*)Y->data;
475ef66eb69SBarry Smith 
476ef66eb69SBarry Smith   PetscFunctionBegin;
4778fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
478ef66eb69SBarry Smith     shell->vshift = 0.0;
479ef66eb69SBarry Smith     shell->vscale = 1.0;
4800298fd71SBarry Smith     shell->dshift = NULL;
4810298fd71SBarry Smith     shell->left   = NULL;
4820298fd71SBarry Smith     shell->right  = NULL;
4837fb7d96aSJed Brown     if (shell->mult) {
484ef66eb69SBarry Smith       Y->ops->mult = shell->mult;
4850298fd71SBarry Smith       shell->mult  = NULL;
4867fb7d96aSJed Brown     }
4877fb7d96aSJed Brown     if (shell->multtranspose) {
48818be62a5SSatish Balay       Y->ops->multtranspose = shell->multtranspose;
4890298fd71SBarry Smith       shell->multtranspose  = NULL;
4907fb7d96aSJed Brown     }
4917fb7d96aSJed Brown     if (shell->getdiagonal) {
49218be62a5SSatish Balay       Y->ops->getdiagonal = shell->getdiagonal;
4930298fd71SBarry Smith       shell->getdiagonal  = NULL;
4947fb7d96aSJed Brown     }
4957fb7d96aSJed Brown     shell->usingscaled = PETSC_FALSE;
496ef66eb69SBarry Smith   }
497ef66eb69SBarry Smith   PetscFunctionReturn(0);
498ef66eb69SBarry Smith }
499ef66eb69SBarry Smith 
500cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
501b951964fSBarry Smith 
5023b49f96aSBarry Smith #undef __FUNCT__
5033b49f96aSBarry Smith #define __FUNCT__ "MatMissingDiagonal_Shell"
5043b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
5053b49f96aSBarry Smith {
5063b49f96aSBarry Smith   PetscFunctionBegin;
5073b49f96aSBarry Smith   *missing = PETSC_FALSE;
5083b49f96aSBarry Smith   PetscFunctionReturn(0);
5093b49f96aSBarry Smith }
5103b49f96aSBarry Smith 
51109dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
51220563c6bSBarry Smith                                        0,
51320563c6bSBarry Smith                                        0,
51420563c6bSBarry Smith                                        0,
51597304618SKris Buschelman                                 /* 4*/ 0,
51620563c6bSBarry Smith                                        0,
517b951964fSBarry Smith                                        0,
518b951964fSBarry Smith                                        0,
519b951964fSBarry Smith                                        0,
520b951964fSBarry Smith                                        0,
52197304618SKris Buschelman                                 /*10*/ 0,
522b951964fSBarry Smith                                        0,
523b951964fSBarry Smith                                        0,
524b951964fSBarry Smith                                        0,
525b951964fSBarry Smith                                        0,
52697304618SKris Buschelman                                 /*15*/ 0,
527b951964fSBarry Smith                                        0,
528b951964fSBarry Smith                                        0,
5298fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
530b951964fSBarry Smith                                        0,
53197304618SKris Buschelman                                 /*20*/ 0,
532ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
533b951964fSBarry Smith                                        0,
534b951964fSBarry Smith                                        0,
535d519adbfSMatthew Knepley                                 /*24*/ 0,
536b951964fSBarry Smith                                        0,
537b951964fSBarry Smith                                        0,
538b951964fSBarry Smith                                        0,
539b951964fSBarry Smith                                        0,
540d519adbfSMatthew Knepley                                 /*29*/ 0,
541b951964fSBarry Smith                                        0,
542273d9f13SBarry Smith                                        0,
543b951964fSBarry Smith                                        0,
544b951964fSBarry Smith                                        0,
545d519adbfSMatthew Knepley                                 /*34*/ 0,
546b951964fSBarry Smith                                        0,
547b951964fSBarry Smith                                        0,
54809dc0095SBarry Smith                                        0,
54909dc0095SBarry Smith                                        0,
550d519adbfSMatthew Knepley                                 /*39*/ 0,
55109dc0095SBarry Smith                                        0,
55209dc0095SBarry Smith                                        0,
55309dc0095SBarry Smith                                        0,
554*7fabda3fSAlex Fikl                                        MatCopy_Shell,
555d519adbfSMatthew Knepley                                 /*44*/ 0,
556ef66eb69SBarry Smith                                        MatScale_Shell,
557ef66eb69SBarry Smith                                        MatShift_Shell,
5586464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
55909dc0095SBarry Smith                                        0,
560f73d5cc4SBarry Smith                                 /*49*/ 0,
56109dc0095SBarry Smith                                        0,
56209dc0095SBarry Smith                                        0,
56309dc0095SBarry Smith                                        0,
56409dc0095SBarry Smith                                        0,
565d519adbfSMatthew Knepley                                 /*54*/ 0,
56609dc0095SBarry Smith                                        0,
56709dc0095SBarry Smith                                        0,
56809dc0095SBarry Smith                                        0,
56909dc0095SBarry Smith                                        0,
570d519adbfSMatthew Knepley                                 /*59*/ 0,
571b9b97703SBarry Smith                                        MatDestroy_Shell,
57209dc0095SBarry Smith                                        0,
573357abbc8SBarry Smith                                        0,
574273d9f13SBarry Smith                                        0,
575d519adbfSMatthew Knepley                                 /*64*/ 0,
576273d9f13SBarry Smith                                        0,
577273d9f13SBarry Smith                                        0,
578273d9f13SBarry Smith                                        0,
579273d9f13SBarry Smith                                        0,
580d519adbfSMatthew Knepley                                 /*69*/ 0,
581273d9f13SBarry Smith                                        0,
582c87e5d42SMatthew Knepley                                        MatConvert_Shell,
583273d9f13SBarry Smith                                        0,
58497304618SKris Buschelman                                        0,
585d519adbfSMatthew Knepley                                 /*74*/ 0,
58697304618SKris Buschelman                                        0,
58797304618SKris Buschelman                                        0,
58897304618SKris Buschelman                                        0,
58997304618SKris Buschelman                                        0,
590d519adbfSMatthew Knepley                                 /*79*/ 0,
59197304618SKris Buschelman                                        0,
59297304618SKris Buschelman                                        0,
59397304618SKris Buschelman                                        0,
594865e5f61SKris Buschelman                                        0,
595d519adbfSMatthew Knepley                                 /*84*/ 0,
596865e5f61SKris Buschelman                                        0,
597865e5f61SKris Buschelman                                        0,
598865e5f61SKris Buschelman                                        0,
599865e5f61SKris Buschelman                                        0,
600d519adbfSMatthew Knepley                                 /*89*/ 0,
601865e5f61SKris Buschelman                                        0,
602865e5f61SKris Buschelman                                        0,
603865e5f61SKris Buschelman                                        0,
604865e5f61SKris Buschelman                                        0,
605d519adbfSMatthew Knepley                                 /*94*/ 0,
606865e5f61SKris Buschelman                                        0,
607865e5f61SKris Buschelman                                        0,
6083964eb88SJed Brown                                        0,
6093964eb88SJed Brown                                        0,
6103964eb88SJed Brown                                 /*99*/ 0,
6113964eb88SJed Brown                                        0,
6123964eb88SJed Brown                                        0,
6133964eb88SJed Brown                                        0,
6143964eb88SJed Brown                                        0,
6153964eb88SJed Brown                                /*104*/ 0,
6163964eb88SJed Brown                                        0,
6173964eb88SJed Brown                                        0,
6183964eb88SJed Brown                                        0,
6193964eb88SJed Brown                                        0,
6203964eb88SJed Brown                                /*109*/ 0,
6213964eb88SJed Brown                                        0,
6223964eb88SJed Brown                                        0,
6233964eb88SJed Brown                                        0,
6243b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
6253964eb88SJed Brown                                /*114*/ 0,
6263964eb88SJed Brown                                        0,
6273964eb88SJed Brown                                        0,
6283964eb88SJed Brown                                        0,
6293964eb88SJed Brown                                        0,
6303964eb88SJed Brown                                /*119*/ 0,
6313964eb88SJed Brown                                        0,
6323964eb88SJed Brown                                        0,
6333964eb88SJed Brown                                        0,
6343964eb88SJed Brown                                        0,
6353964eb88SJed Brown                                /*124*/ 0,
6363964eb88SJed Brown                                        0,
6373964eb88SJed Brown                                        0,
6383964eb88SJed Brown                                        0,
6393964eb88SJed Brown                                        0,
6403964eb88SJed Brown                                /*129*/ 0,
6413964eb88SJed Brown                                        0,
6423964eb88SJed Brown                                        0,
6433964eb88SJed Brown                                        0,
6443964eb88SJed Brown                                        0,
6453964eb88SJed Brown                                /*134*/ 0,
6463964eb88SJed Brown                                        0,
6473964eb88SJed Brown                                        0,
6483964eb88SJed Brown                                        0,
6493964eb88SJed Brown                                        0,
6503964eb88SJed Brown                                /*139*/ 0,
651f9426fe0SMark Adams                                        0,
6523964eb88SJed Brown                                        0
6533964eb88SJed Brown };
654273d9f13SBarry Smith 
6550bad9183SKris Buschelman /*MC
656fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
6570bad9183SKris Buschelman 
6580bad9183SKris Buschelman   Level: advanced
6590bad9183SKris Buschelman 
6600bad9183SKris Buschelman .seealso: MatCreateShell
6610bad9183SKris Buschelman M*/
6620bad9183SKris Buschelman 
663711e205bSSatish Balay #undef __FUNCT__
664711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell"
6658cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
666273d9f13SBarry Smith {
667273d9f13SBarry Smith   Mat_Shell      *b;
668dfbe8321SBarry Smith   PetscErrorCode ierr;
669273d9f13SBarry Smith 
670273d9f13SBarry Smith   PetscFunctionBegin;
671273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
672273d9f13SBarry Smith 
673b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
674273d9f13SBarry Smith   A->data = (void*)b;
675273d9f13SBarry Smith 
67626283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
67726283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
678273d9f13SBarry Smith 
679273d9f13SBarry Smith   b->ctx           = 0;
680ef66eb69SBarry Smith   b->vshift        = 0.0;
681ef66eb69SBarry Smith   b->vscale        = 1.0;
682ef66eb69SBarry Smith   b->mult          = 0;
68318be62a5SSatish Balay   b->multtranspose = 0;
68418be62a5SSatish Balay   b->getdiagonal   = 0;
685273d9f13SBarry Smith   A->assembled     = PETSC_TRUE;
686210f0121SBarry Smith   A->preallocated  = PETSC_FALSE;
6872205254eSKarl Rupp 
68817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
689273d9f13SBarry Smith   PetscFunctionReturn(0);
690273d9f13SBarry Smith }
691e51e0e81SBarry Smith 
692711e205bSSatish Balay #undef __FUNCT__
693711e205bSSatish Balay #define __FUNCT__ "MatCreateShell"
6944b828684SBarry Smith /*@C
695052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
696ff756334SLois Curfman McInnes    private data storage format.
697e51e0e81SBarry Smith 
698c7fcc2eaSBarry Smith   Collective on MPI_Comm
699c7fcc2eaSBarry Smith 
700e51e0e81SBarry Smith    Input Parameters:
701c7fcc2eaSBarry Smith +  comm - MPI communicator
702c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
703c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
704c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
705c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
706c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
707e51e0e81SBarry Smith 
708ff756334SLois Curfman McInnes    Output Parameter:
70944cd7ae7SLois Curfman McInnes .  A - the matrix
710e51e0e81SBarry Smith 
711ff2fd236SBarry Smith    Level: advanced
712ff2fd236SBarry Smith 
713f39d1f56SLois Curfman McInnes   Usage:
7147b2a1423SBarry Smith $    extern int mult(Mat,Vec,Vec);
715f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
716c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
717f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
718f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
719f39d1f56SLois Curfman McInnes 
720ff756334SLois Curfman McInnes    Notes:
721ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
722ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
723ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
724e51e0e81SBarry Smith 
725daf670e6SBarry Smith    Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this
726daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
727daf670e6SBarry Smith     in as the ctx argument.
728f38a8d56SBarry Smith 
729f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
730f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
731645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
732645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
733645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
734645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
735645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
736645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
737645985a0SLois Curfman McInnes    For example,
738f39d1f56SLois Curfman McInnes 
739f39d1f56SLois Curfman McInnes $
740f39d1f56SLois Curfman McInnes $     Vec x, y
7417b2a1423SBarry Smith $     extern int mult(Mat,Vec,Vec);
742645985a0SLois Curfman McInnes $     Mat A
743f39d1f56SLois Curfman McInnes $
744c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
745c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
746f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
747c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
748c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
749c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
750645985a0SLois Curfman McInnes $     MatMult(A,x,y);
751645985a0SLois Curfman McInnes $     MatDestroy(A);
752f39d1f56SLois Curfman McInnes $     VecDestroy(y); VecDestroy(x);
753645985a0SLois Curfman McInnes $
754e51e0e81SBarry Smith 
7550b627109SLois Curfman McInnes .keywords: matrix, shell, create
7560b627109SLois Curfman McInnes 
757ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
758e51e0e81SBarry Smith @*/
7597087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
760e51e0e81SBarry Smith {
761dfbe8321SBarry Smith   PetscErrorCode ierr;
762ed3cc1f0SBarry Smith 
7633a40ed3dSBarry Smith   PetscFunctionBegin;
764f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
765f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
766273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
767273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
7680fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
769273d9f13SBarry Smith   PetscFunctionReturn(0);
770c7fcc2eaSBarry Smith }
771c7fcc2eaSBarry Smith 
772711e205bSSatish Balay #undef __FUNCT__
773711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext"
774c6866cfdSSatish Balay /*@
775273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
776c7fcc2eaSBarry Smith 
7773f9fe445SBarry Smith    Logically Collective on Mat
778c7fcc2eaSBarry Smith 
779273d9f13SBarry Smith     Input Parameters:
780273d9f13SBarry Smith +   mat - the shell matrix
781273d9f13SBarry Smith -   ctx - the context
782273d9f13SBarry Smith 
783273d9f13SBarry Smith    Level: advanced
784273d9f13SBarry Smith 
785daf670e6SBarry Smith    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
786daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
787273d9f13SBarry Smith 
788273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
7890bc0a719SBarry Smith @*/
7907087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
791273d9f13SBarry Smith {
792273d9f13SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
793dfbe8321SBarry Smith   PetscErrorCode ierr;
794ace3abfcSBarry Smith   PetscBool      flg;
795273d9f13SBarry Smith 
796273d9f13SBarry Smith   PetscFunctionBegin;
7970700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
798251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
799273d9f13SBarry Smith   if (flg) {
800273d9f13SBarry Smith     shell->ctx = ctx;
801ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
8023a40ed3dSBarry Smith   PetscFunctionReturn(0);
803e51e0e81SBarry Smith }
804e51e0e81SBarry Smith 
805711e205bSSatish Balay #undef __FUNCT__
806711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation"
807c16cb8f2SBarry Smith /*@C
8083a3eedf2SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
8093a3eedf2SBarry Smith                            a shell matrix.
810e51e0e81SBarry Smith 
8113f9fe445SBarry Smith    Logically Collective on Mat
812fee21e36SBarry Smith 
813c7fcc2eaSBarry Smith     Input Parameters:
814c7fcc2eaSBarry Smith +   mat - the shell matrix
815c7fcc2eaSBarry Smith .   op - the name of the operation
816c7fcc2eaSBarry Smith -   f - the function that provides the operation.
817c7fcc2eaSBarry Smith 
81815091d37SBarry Smith    Level: advanced
81915091d37SBarry Smith 
820fae171e0SBarry Smith     Usage:
821e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
822f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
823c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
8240b627109SLois Curfman McInnes 
825a62d957aSLois Curfman McInnes     Notes:
826e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
8271c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
828a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
8291c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
830a62d957aSLois Curfman McInnes 
83125296bd5SBarry Smith     All user-provided functions (execept for MATOP_DESTROY) should have the same calling
832deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
833deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
834deebb3c3SLois Curfman McInnes     routines, e.g.,
835a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
836a62d957aSLois Curfman McInnes 
8374aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
8384aa34b0aSBarry Smith     nonzero on failure.
8394aa34b0aSBarry Smith 
840a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
841a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
842a62d957aSLois Curfman McInnes     set by MatCreateShell().
843a62d957aSLois Curfman McInnes 
8442a7a6963SBarry Smith     Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
845501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
846501d9185SBarry Smith 
847a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation
848a62d957aSLois Curfman McInnes 
849ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
850e51e0e81SBarry Smith @*/
8517087cfbeSBarry Smith PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
852e51e0e81SBarry Smith {
853dfbe8321SBarry Smith   PetscErrorCode ierr;
854ace3abfcSBarry Smith   PetscBool      flg;
855273d9f13SBarry Smith 
8563a40ed3dSBarry Smith   PetscFunctionBegin;
8570700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
8585edf6516SJed Brown   switch (op) {
8595edf6516SJed Brown   case MATOP_DESTROY:
860251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
861273d9f13SBarry Smith     if (flg) {
862a62d957aSLois Curfman McInnes       Mat_Shell *shell = (Mat_Shell*)mat->data;
8636849ba73SBarry Smith       shell->destroy = (PetscErrorCode (*)(Mat))f;
8646849ba73SBarry Smith     } else mat->ops->destroy = (PetscErrorCode (*)(Mat))f;
8655edf6516SJed Brown     break;
8666464bf51SAlex Fikl   case MATOP_DIAGONAL_SET:
8676464bf51SAlex Fikl     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
8686464bf51SAlex Fikl     if (flg) {
8696464bf51SAlex Fikl       Mat_Shell *shell = (Mat_Shell*)mat->data;
8706464bf51SAlex Fikl       shell->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f;
8716464bf51SAlex Fikl     } else {
8726464bf51SAlex Fikl       mat->ops->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f;
8736464bf51SAlex Fikl     }
8746464bf51SAlex Fikl     break;
8755edf6516SJed Brown   case MATOP_VIEW:
8765edf6516SJed Brown     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
8775edf6516SJed Brown     break;
8785edf6516SJed Brown   case MATOP_MULT:
8795edf6516SJed Brown     mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
8805edf6516SJed Brown     if (!mat->ops->multadd) mat->ops->multadd = MatMultAdd_Shell;
8815edf6516SJed Brown     break;
8825edf6516SJed Brown   case MATOP_MULT_TRANSPOSE:
8835edf6516SJed Brown     mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
8845edf6516SJed Brown     if (!mat->ops->multtransposeadd) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell;
8855edf6516SJed Brown     break;
8865edf6516SJed Brown   default:
8875edf6516SJed Brown     (((void(**)(void))mat->ops)[op]) = f;
888a62d957aSLois Curfman McInnes   }
8893a40ed3dSBarry Smith   PetscFunctionReturn(0);
890e51e0e81SBarry Smith }
891f0479e8cSBarry Smith 
892711e205bSSatish Balay #undef __FUNCT__
893711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation"
894d4bb536fSBarry Smith /*@C
895d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
896d4bb536fSBarry Smith 
897c7fcc2eaSBarry Smith     Not Collective
898c7fcc2eaSBarry Smith 
899d4bb536fSBarry Smith     Input Parameters:
900c7fcc2eaSBarry Smith +   mat - the shell matrix
901c7fcc2eaSBarry Smith -   op - the name of the operation
902d4bb536fSBarry Smith 
903d4bb536fSBarry Smith     Output Parameter:
904d4bb536fSBarry Smith .   f - the function that provides the operation.
905d4bb536fSBarry Smith 
90615091d37SBarry Smith     Level: advanced
90715091d37SBarry Smith 
908d4bb536fSBarry Smith     Notes:
909e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
910d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
911d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
912d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
913d4bb536fSBarry Smith 
914d4bb536fSBarry Smith     All user-provided functions have the same calling
915d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
916d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
917d4bb536fSBarry Smith     routines, e.g.,
918d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
919d4bb536fSBarry Smith 
920d4bb536fSBarry Smith     Within each user-defined routine, the user should call
921d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
922d4bb536fSBarry Smith     set by MatCreateShell().
923d4bb536fSBarry Smith 
924d4bb536fSBarry Smith .keywords: matrix, shell, set, operation
925d4bb536fSBarry Smith 
926ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
927d4bb536fSBarry Smith @*/
9287087cfbeSBarry Smith PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
929d4bb536fSBarry Smith {
930dfbe8321SBarry Smith   PetscErrorCode ierr;
931ace3abfcSBarry Smith   PetscBool      flg;
932273d9f13SBarry Smith 
9333a40ed3dSBarry Smith   PetscFunctionBegin;
9340700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
935d4bb536fSBarry Smith   if (op == MATOP_DESTROY) {
936251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
937273d9f13SBarry Smith     if (flg) {
938d4bb536fSBarry Smith       Mat_Shell *shell = (Mat_Shell*)mat->data;
939c134de8dSSatish Balay       *f = (void (*)(void))shell->destroy;
940c7fcc2eaSBarry Smith     } else {
941c134de8dSSatish Balay       *f = (void (*)(void))mat->ops->destroy;
942d4bb536fSBarry Smith     }
943c7fcc2eaSBarry Smith   } else if (op == MATOP_VIEW) {
944c134de8dSSatish Balay     *f = (void (*)(void))mat->ops->view;
945c7fcc2eaSBarry Smith   } else {
946c134de8dSSatish Balay     *f = (((void (**)(void))mat->ops)[op]);
947d4bb536fSBarry Smith   }
9483a40ed3dSBarry Smith   PetscFunctionReturn(0);
949d4bb536fSBarry Smith }
950d4bb536fSBarry Smith 
951