xref: /petsc/src/mat/impls/shell/shell.c (revision ef5c7bd2f4ec6f9849104a8422b5777b3b2f72fb)
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*/
945960306SStefano Zampini #include <petsc/private/vecimpl.h>
10e51e0e81SBarry Smith 
1128f357e3SAlex Fikl struct _MatShellOps {
12976e8c5aSLisandro Dalcin   /*  3 */ PetscErrorCode (*mult)(Mat,Vec,Vec);
13976e8c5aSLisandro Dalcin   /*  5 */ PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
14976e8c5aSLisandro Dalcin   /* 17 */ PetscErrorCode (*getdiagonal)(Mat,Vec);
15976e8c5aSLisandro Dalcin   /* 43 */ PetscErrorCode (*copy)(Mat,Mat,MatStructure);
16976e8c5aSLisandro Dalcin   /* 60 */ PetscErrorCode (*destroy)(Mat);
1728f357e3SAlex Fikl };
1828f357e3SAlex Fikl 
1928f357e3SAlex Fikl typedef struct {
2028f357e3SAlex Fikl   struct _MatShellOps ops[1];
212205254eSKarl Rupp 
22ef66eb69SBarry Smith   PetscScalar vscale,vshift;
238fe8eb27SJed Brown   Vec         dshift;
248fe8eb27SJed Brown   Vec         left,right;
258fe8eb27SJed Brown   Vec         left_work,right_work;
265edf6516SJed Brown   Vec         left_add_work,right_add_work;
27*ef5c7bd2SStefano Zampini   /* support for MatAXPY */
289f137db4SBarry Smith   Mat              axpy;
299f137db4SBarry Smith   PetscScalar      axpy_vscale;
30*ef5c7bd2SStefano Zampini   PetscObjectState axpy_state;
310c0fd78eSBarry Smith   PetscBool        managescalingshifts; /* The user will manage the scaling and shifts for the MATSHELL, not the default */
3245960306SStefano Zampini   /* support for ZeroRows/Columns operations */
3345960306SStefano Zampini   IS          zrows;
3445960306SStefano Zampini   IS          zcols;
3545960306SStefano Zampini   Vec         zvals;
3645960306SStefano Zampini   Vec         zvals_w;
3745960306SStefano Zampini   VecScatter  zvals_sct_r;
3845960306SStefano Zampini   VecScatter  zvals_sct_c;
3920563c6bSBarry Smith   void        *ctx;
4088cf3e7dSBarry Smith } Mat_Shell;
410c0fd78eSBarry Smith 
4245960306SStefano Zampini 
4345960306SStefano Zampini /*
4445960306SStefano Zampini      Store and scale values on zeroed rows
4545960306SStefano Zampini      xx = [x_1, 0], 0 on zeroed columns
4645960306SStefano Zampini */
4745960306SStefano Zampini static PetscErrorCode MatShellPreZeroRight(Mat A,Vec x,Vec *xx)
4845960306SStefano Zampini {
4945960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
5045960306SStefano Zampini   PetscErrorCode ierr;
5145960306SStefano Zampini 
5245960306SStefano Zampini   PetscFunctionBegin;
5345960306SStefano Zampini   *xx = x;
5445960306SStefano Zampini   if (shell->zrows) {
5545960306SStefano Zampini     ierr = VecSet(shell->zvals_w,0.0);CHKERRQ(ierr);
5645960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5745960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5845960306SStefano Zampini     ierr = VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);CHKERRQ(ierr);
5945960306SStefano Zampini   }
6045960306SStefano Zampini   if (shell->zcols) {
6145960306SStefano Zampini     if (!shell->right_work) {
6245960306SStefano Zampini       ierr = MatCreateVecs(A,&shell->right_work,NULL);CHKERRQ(ierr);
6345960306SStefano Zampini     }
6445960306SStefano Zampini     ierr = VecCopy(x,shell->right_work);CHKERRQ(ierr);
6545960306SStefano Zampini     ierr = VecISSet(shell->right_work,shell->zcols,0.0);CHKERRQ(ierr);
6645960306SStefano Zampini     *xx  = shell->right_work;
6745960306SStefano Zampini   }
6845960306SStefano Zampini   PetscFunctionReturn(0);
6945960306SStefano Zampini }
7045960306SStefano Zampini 
7145960306SStefano Zampini /* Insert properly diagonally scaled values stored in MatShellPreZeroRight */
7245960306SStefano Zampini static PetscErrorCode MatShellPostZeroLeft(Mat A,Vec x)
7345960306SStefano Zampini {
7445960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
7545960306SStefano Zampini   PetscErrorCode ierr;
7645960306SStefano Zampini 
7745960306SStefano Zampini   PetscFunctionBegin;
7845960306SStefano Zampini   if (shell->zrows) {
7945960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8045960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8145960306SStefano Zampini   }
8245960306SStefano Zampini   PetscFunctionReturn(0);
8345960306SStefano Zampini }
8445960306SStefano Zampini 
8545960306SStefano Zampini /*
8645960306SStefano Zampini      Store and scale values on zeroed rows
8745960306SStefano Zampini      xx = [x_1, 0], 0 on zeroed rows
8845960306SStefano Zampini */
8945960306SStefano Zampini static PetscErrorCode MatShellPreZeroLeft(Mat A,Vec x,Vec *xx)
9045960306SStefano Zampini {
9145960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
9245960306SStefano Zampini   PetscErrorCode ierr;
9345960306SStefano Zampini 
9445960306SStefano Zampini   PetscFunctionBegin;
9545960306SStefano Zampini   *xx = NULL;
9645960306SStefano Zampini   if (!shell->zrows) {
9745960306SStefano Zampini     *xx = x;
9845960306SStefano Zampini   } else {
9945960306SStefano Zampini     if (!shell->left_work) {
10045960306SStefano Zampini       ierr = MatCreateVecs(A,NULL,&shell->left_work);CHKERRQ(ierr);
10145960306SStefano Zampini     }
10245960306SStefano Zampini     ierr = VecCopy(x,shell->left_work);CHKERRQ(ierr);
10345960306SStefano Zampini     ierr = VecSet(shell->zvals_w,0.0);CHKERRQ(ierr);
10445960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10545960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10645960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10745960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10845960306SStefano Zampini     ierr = VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);CHKERRQ(ierr);
10945960306SStefano Zampini     *xx  = shell->left_work;
11045960306SStefano Zampini   }
11145960306SStefano Zampini   PetscFunctionReturn(0);
11245960306SStefano Zampini }
11345960306SStefano Zampini 
11445960306SStefano Zampini /* Zero zero-columns contributions, sum contributions from properly scaled values stored in MatShellPreZeroLeft */
11545960306SStefano Zampini static PetscErrorCode MatShellPostZeroRight(Mat A,Vec x)
11645960306SStefano Zampini {
11745960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
11845960306SStefano Zampini   PetscErrorCode ierr;
11945960306SStefano Zampini 
12045960306SStefano Zampini   PetscFunctionBegin;
12145960306SStefano Zampini   if (shell->zcols) {
12245960306SStefano Zampini     ierr = VecISSet(x,shell->zcols,0.0);CHKERRQ(ierr);
12345960306SStefano Zampini   }
12445960306SStefano Zampini   if (shell->zrows) {
12545960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12645960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12745960306SStefano Zampini   }
12845960306SStefano Zampini   PetscFunctionReturn(0);
12945960306SStefano Zampini }
13045960306SStefano Zampini 
1318fe8eb27SJed Brown /*
1320c0fd78eSBarry Smith       xx = diag(left)*x
1338fe8eb27SJed Brown */
1348fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
1358fe8eb27SJed Brown {
1368fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1378fe8eb27SJed Brown   PetscErrorCode ierr;
1388fe8eb27SJed Brown 
1398fe8eb27SJed Brown   PetscFunctionBegin;
1400298fd71SBarry Smith   *xx = NULL;
1418fe8eb27SJed Brown   if (!shell->left) {
1428fe8eb27SJed Brown     *xx = x;
1438fe8eb27SJed Brown   } else {
1448fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
1458fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
1468fe8eb27SJed Brown     *xx  = shell->left_work;
1478fe8eb27SJed Brown   }
1488fe8eb27SJed Brown   PetscFunctionReturn(0);
1498fe8eb27SJed Brown }
1508fe8eb27SJed Brown 
1510c0fd78eSBarry Smith /*
1520c0fd78eSBarry Smith      xx = diag(right)*x
1530c0fd78eSBarry Smith */
1548fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
1558fe8eb27SJed Brown {
1568fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1578fe8eb27SJed Brown   PetscErrorCode ierr;
1588fe8eb27SJed Brown 
1598fe8eb27SJed Brown   PetscFunctionBegin;
1600298fd71SBarry Smith   *xx = NULL;
1618fe8eb27SJed Brown   if (!shell->right) {
1628fe8eb27SJed Brown     *xx = x;
1638fe8eb27SJed Brown   } else {
1648fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
1658fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
1668fe8eb27SJed Brown     *xx  = shell->right_work;
1678fe8eb27SJed Brown   }
1688fe8eb27SJed Brown   PetscFunctionReturn(0);
1698fe8eb27SJed Brown }
1708fe8eb27SJed Brown 
1710c0fd78eSBarry Smith /*
1720c0fd78eSBarry Smith     x = diag(left)*x
1730c0fd78eSBarry Smith */
1748fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
1758fe8eb27SJed Brown {
1768fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1778fe8eb27SJed Brown   PetscErrorCode ierr;
1788fe8eb27SJed Brown 
1798fe8eb27SJed Brown   PetscFunctionBegin;
1808fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
1818fe8eb27SJed Brown   PetscFunctionReturn(0);
1828fe8eb27SJed Brown }
1838fe8eb27SJed Brown 
1840c0fd78eSBarry Smith /*
1850c0fd78eSBarry Smith     x = diag(right)*x
1860c0fd78eSBarry Smith */
1878fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
1888fe8eb27SJed Brown {
1898fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1908fe8eb27SJed Brown   PetscErrorCode ierr;
1918fe8eb27SJed Brown 
1928fe8eb27SJed Brown   PetscFunctionBegin;
1938fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
1948fe8eb27SJed Brown   PetscFunctionReturn(0);
1958fe8eb27SJed Brown }
1968fe8eb27SJed Brown 
1970c0fd78eSBarry Smith /*
1980c0fd78eSBarry Smith          Y = vscale*Y + diag(dshift)*X + vshift*X
1990c0fd78eSBarry Smith 
2000c0fd78eSBarry Smith          On input Y already contains A*x
2010c0fd78eSBarry Smith */
2028fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
2038fe8eb27SJed Brown {
2048fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
2058fe8eb27SJed Brown   PetscErrorCode ierr;
2068fe8eb27SJed Brown 
2078fe8eb27SJed Brown   PetscFunctionBegin;
2088fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
2098fe8eb27SJed Brown     PetscInt          i,m;
2108fe8eb27SJed Brown     const PetscScalar *x,*d;
2118fe8eb27SJed Brown     PetscScalar       *y;
2128fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
2138fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
2148fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
2158fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
2168fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
2178fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
2188fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
2198fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
2200c0fd78eSBarry Smith   } else {
221027c5fb5SJed Brown     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
2228fe8eb27SJed Brown   }
223d4c7de66SBarry Smith   if (shell->vshift != 0.0) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */
2248fe8eb27SJed Brown   PetscFunctionReturn(0);
2258fe8eb27SJed Brown }
226e51e0e81SBarry Smith 
227789d8953SBarry Smith PetscErrorCode MatShellGetContext_Shell(Mat mat,void *ctx)
228789d8953SBarry Smith {
229789d8953SBarry Smith   PetscFunctionBegin;
230789d8953SBarry Smith   *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
231789d8953SBarry Smith   PetscFunctionReturn(0);
232789d8953SBarry Smith }
233789d8953SBarry Smith 
2349d225801SJed Brown /*@
235a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
236b4fd4287SBarry Smith 
237c7fcc2eaSBarry Smith     Not Collective
238c7fcc2eaSBarry Smith 
239b4fd4287SBarry Smith     Input Parameter:
240b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
241b4fd4287SBarry Smith 
242b4fd4287SBarry Smith     Output Parameter:
243b4fd4287SBarry Smith .   ctx - the user provided context
244b4fd4287SBarry Smith 
24515091d37SBarry Smith     Level: advanced
24615091d37SBarry Smith 
24795452b02SPatrick Sanan    Fortran Notes:
24895452b02SPatrick Sanan     To use this from Fortran you must write a Fortran interface definition for this
249daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
250a62d957aSLois Curfman McInnes 
251ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
2520bc0a719SBarry Smith @*/
2538fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
254b4fd4287SBarry Smith {
255dfbe8321SBarry Smith   PetscErrorCode ierr;
256273d9f13SBarry Smith 
2573a40ed3dSBarry Smith   PetscFunctionBegin;
2580700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2594482741eSBarry Smith   PetscValidPointer(ctx,2);
260789d8953SBarry Smith   ierr = PetscUseMethod(mat,"MatShellGetContext_C",(Mat,void*),(mat,ctx));CHKERRQ(ierr);
2613a40ed3dSBarry Smith   PetscFunctionReturn(0);
262b4fd4287SBarry Smith }
263b4fd4287SBarry Smith 
26445960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Local_Shell(Mat mat,PetscInt nr,PetscInt rows[],PetscInt nc,PetscInt cols[],PetscScalar diag,PetscBool rc)
26545960306SStefano Zampini {
26645960306SStefano Zampini   PetscErrorCode ierr;
26745960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)mat->data;
26845960306SStefano Zampini   Vec            x = NULL,b = NULL;
26945960306SStefano Zampini   IS             is1, is2;
27045960306SStefano Zampini   const PetscInt *ridxs;
27145960306SStefano Zampini   PetscInt       *idxs,*gidxs;
27245960306SStefano Zampini   PetscInt       cum,rst,cst,i;
27345960306SStefano Zampini 
27445960306SStefano Zampini   PetscFunctionBegin;
27545960306SStefano Zampini   if (!shell->zvals) {
27645960306SStefano Zampini     ierr = MatCreateVecs(mat,NULL,&shell->zvals);CHKERRQ(ierr);
27745960306SStefano Zampini   }
27845960306SStefano Zampini   if (!shell->zvals_w) {
27945960306SStefano Zampini     ierr = VecDuplicate(shell->zvals,&shell->zvals_w);CHKERRQ(ierr);
28045960306SStefano Zampini   }
28145960306SStefano Zampini   ierr = MatGetOwnershipRange(mat,&rst,NULL);CHKERRQ(ierr);
28245960306SStefano Zampini   ierr = MatGetOwnershipRangeColumn(mat,&cst,NULL);CHKERRQ(ierr);
28345960306SStefano Zampini 
28445960306SStefano Zampini   /* Expand/create index set of zeroed rows */
28545960306SStefano Zampini   ierr = PetscMalloc1(nr,&idxs);CHKERRQ(ierr);
28645960306SStefano Zampini   for (i = 0; i < nr; i++) idxs[i] = rows[i] + rst;
28745960306SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,nr,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr);
28845960306SStefano Zampini   ierr = ISSort(is1);CHKERRQ(ierr);
28945960306SStefano Zampini   ierr = VecISSet(shell->zvals,is1,diag);CHKERRQ(ierr);
29045960306SStefano Zampini   if (shell->zrows) {
29145960306SStefano Zampini     ierr = ISSum(shell->zrows,is1,&is2);CHKERRQ(ierr);
29245960306SStefano Zampini     ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr);
29345960306SStefano Zampini     ierr = ISDestroy(&is1);CHKERRQ(ierr);
29445960306SStefano Zampini     shell->zrows = is2;
29545960306SStefano Zampini   } else shell->zrows = is1;
29645960306SStefano Zampini 
29745960306SStefano Zampini   /* Create scatters for diagonal values communications */
29845960306SStefano Zampini   ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr);
29945960306SStefano Zampini   ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr);
30045960306SStefano Zampini 
30145960306SStefano Zampini   /* row scatter: from/to left vector */
30245960306SStefano Zampini   ierr = MatCreateVecs(mat,&x,&b);CHKERRQ(ierr);
30345960306SStefano Zampini   ierr = VecScatterCreate(b,shell->zrows,shell->zvals_w,shell->zrows,&shell->zvals_sct_r);CHKERRQ(ierr);
30445960306SStefano Zampini 
30545960306SStefano Zampini   /* col scatter: from right vector to left vector */
30645960306SStefano Zampini   ierr = ISGetIndices(shell->zrows,&ridxs);CHKERRQ(ierr);
30745960306SStefano Zampini   ierr = ISGetLocalSize(shell->zrows,&nr);CHKERRQ(ierr);
30845960306SStefano Zampini   ierr = PetscMalloc1(nr,&gidxs);CHKERRQ(ierr);
30945960306SStefano Zampini   for (i = 0, cum  = 0; i < nr; i++) {
31045960306SStefano Zampini     if (ridxs[i] >= mat->cmap->N) continue;
31145960306SStefano Zampini     gidxs[cum] = ridxs[i];
31245960306SStefano Zampini     cum++;
31345960306SStefano Zampini   }
31445960306SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)mat),cum,gidxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr);
31545960306SStefano Zampini   ierr = VecScatterCreate(x,is1,shell->zvals_w,is1,&shell->zvals_sct_c);CHKERRQ(ierr);
31645960306SStefano Zampini   ierr = ISDestroy(&is1);CHKERRQ(ierr);
31745960306SStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
31845960306SStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
31945960306SStefano Zampini 
32045960306SStefano Zampini   /* Expand/create index set of zeroed columns */
32145960306SStefano Zampini   if (rc) {
32245960306SStefano Zampini     ierr = PetscMalloc1(nc,&idxs);CHKERRQ(ierr);
32345960306SStefano Zampini     for (i = 0; i < nc; i++) idxs[i] = cols[i] + cst;
32445960306SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nc,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr);
32545960306SStefano Zampini     ierr = ISSort(is1);CHKERRQ(ierr);
32645960306SStefano Zampini     if (shell->zcols) {
32745960306SStefano Zampini       ierr = ISSum(shell->zcols,is1,&is2);CHKERRQ(ierr);
32845960306SStefano Zampini       ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr);
32945960306SStefano Zampini       ierr = ISDestroy(&is1);CHKERRQ(ierr);
33045960306SStefano Zampini       shell->zcols = is2;
33145960306SStefano Zampini     } else shell->zcols = is1;
33245960306SStefano Zampini   }
33345960306SStefano Zampini   PetscFunctionReturn(0);
33445960306SStefano Zampini }
33545960306SStefano Zampini 
33645960306SStefano Zampini static PetscErrorCode MatZeroRows_Shell(Mat mat,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
33745960306SStefano Zampini {
338*ef5c7bd2SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)mat->data;
33945960306SStefano Zampini   PetscInt       nr, *lrows;
34045960306SStefano Zampini   PetscErrorCode ierr;
34145960306SStefano Zampini 
34245960306SStefano Zampini   PetscFunctionBegin;
34345960306SStefano Zampini   if (x && b) {
34445960306SStefano Zampini     Vec          xt;
34545960306SStefano Zampini     PetscScalar *vals;
34645960306SStefano Zampini     PetscInt    *gcols,i,st,nl,nc;
34745960306SStefano Zampini 
34845960306SStefano Zampini     ierr = PetscMalloc1(n,&gcols);CHKERRQ(ierr);
34945960306SStefano Zampini     for (i = 0, nc = 0; i < n; i++) if (rows[i] < mat->cmap->N) gcols[nc++] = rows[i];
35045960306SStefano Zampini 
35145960306SStefano Zampini     ierr = MatCreateVecs(mat,&xt,NULL);CHKERRQ(ierr);
35245960306SStefano Zampini     ierr = VecCopy(x,xt);CHKERRQ(ierr);
35345960306SStefano Zampini     ierr = PetscCalloc1(nc,&vals);CHKERRQ(ierr);
35445960306SStefano Zampini     ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */
35545960306SStefano Zampini     ierr = PetscFree(vals);CHKERRQ(ierr);
35645960306SStefano Zampini     ierr = VecAssemblyBegin(xt);CHKERRQ(ierr);
35745960306SStefano Zampini     ierr = VecAssemblyEnd(xt);CHKERRQ(ierr);
35845960306SStefano Zampini     ierr = VecAYPX(xt,-1.0,x);CHKERRQ(ierr);                           /* xt = [0, x2] */
35945960306SStefano Zampini 
36045960306SStefano Zampini     ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr);
36145960306SStefano Zampini     ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr);
36245960306SStefano Zampini     ierr = VecGetArray(xt,&vals);CHKERRQ(ierr);
36345960306SStefano Zampini     for (i = 0; i < nl; i++) {
36445960306SStefano Zampini       PetscInt g = i + st;
36545960306SStefano Zampini       if (g > mat->rmap->N) continue;
36645960306SStefano Zampini       if (PetscAbsScalar(vals[i]) == 0.0) continue;
36745960306SStefano Zampini       ierr = VecSetValue(b,g,diag*vals[i],INSERT_VALUES);CHKERRQ(ierr);
36845960306SStefano Zampini     }
36945960306SStefano Zampini     ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr);
37045960306SStefano Zampini     ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
37145960306SStefano Zampini     ierr = VecAssemblyEnd(b);CHKERRQ(ierr);                            /* b  = [b1, x2 * diag] */
37245960306SStefano Zampini     ierr = VecDestroy(&xt);CHKERRQ(ierr);
37345960306SStefano Zampini     ierr = PetscFree(gcols);CHKERRQ(ierr);
37445960306SStefano Zampini   }
37545960306SStefano Zampini   ierr = PetscLayoutMapLocal(mat->rmap,n,rows,&nr,&lrows,NULL);CHKERRQ(ierr);
37645960306SStefano Zampini   ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,0,NULL,diag,PETSC_FALSE);CHKERRQ(ierr);
377*ef5c7bd2SStefano Zampini   if (shell->axpy) {
378*ef5c7bd2SStefano Zampini     ierr = MatZeroRows(shell->axpy,n,rows,0.0,NULL,NULL);CHKERRQ(ierr);
379*ef5c7bd2SStefano Zampini   }
38045960306SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
38145960306SStefano Zampini   PetscFunctionReturn(0);
38245960306SStefano Zampini }
38345960306SStefano Zampini 
38445960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Shell(Mat mat,PetscInt n,const PetscInt rowscols[],PetscScalar diag,Vec x,Vec b)
38545960306SStefano Zampini {
386*ef5c7bd2SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)mat->data;
38745960306SStefano Zampini   PetscInt       *lrows, *lcols;
38845960306SStefano Zampini   PetscInt       nr, nc;
38945960306SStefano Zampini   PetscBool      congruent;
39045960306SStefano Zampini   PetscErrorCode ierr;
39145960306SStefano Zampini 
39245960306SStefano Zampini   PetscFunctionBegin;
39345960306SStefano Zampini   if (x && b) {
39445960306SStefano Zampini     Vec          xt, bt;
39545960306SStefano Zampini     PetscScalar *vals;
39645960306SStefano Zampini     PetscInt    *grows,*gcols,i,st,nl;
39745960306SStefano Zampini 
39845960306SStefano Zampini     ierr = PetscMalloc2(n,&grows,n,&gcols);CHKERRQ(ierr);
39945960306SStefano Zampini     for (i = 0, nr = 0; i < n; i++) if (rowscols[i] < mat->rmap->N) grows[nr++] = rowscols[i];
40045960306SStefano Zampini     for (i = 0, nc = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) gcols[nc++] = rowscols[i];
40145960306SStefano Zampini     ierr = PetscCalloc1(n,&vals);CHKERRQ(ierr);
40245960306SStefano Zampini 
40345960306SStefano Zampini     ierr = MatCreateVecs(mat,&xt,&bt);CHKERRQ(ierr);
40445960306SStefano Zampini     ierr = VecCopy(x,xt);CHKERRQ(ierr);
40545960306SStefano Zampini     ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */
40645960306SStefano Zampini     ierr = VecAssemblyBegin(xt);CHKERRQ(ierr);
40745960306SStefano Zampini     ierr = VecAssemblyEnd(xt);CHKERRQ(ierr);
40845960306SStefano Zampini     ierr = VecAXPY(xt,-1.0,x);CHKERRQ(ierr);                           /* xt = [0, -x2] */
40945960306SStefano Zampini     ierr = MatMult(mat,xt,bt);CHKERRQ(ierr);                           /* bt = [-A12*x2,-A22*x2] */
41045960306SStefano Zampini     ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* bt = [-A12*x2,0] */
41145960306SStefano Zampini     ierr = VecAssemblyBegin(bt);CHKERRQ(ierr);
41245960306SStefano Zampini     ierr = VecAssemblyEnd(bt);CHKERRQ(ierr);
41345960306SStefano Zampini     ierr = VecAXPY(b,1.0,bt);CHKERRQ(ierr);                            /* b  = [b1 - A12*x2, b2] */
41445960306SStefano Zampini     ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* b  = [b1 - A12*x2, 0] */
41545960306SStefano Zampini     ierr = VecAssemblyBegin(bt);CHKERRQ(ierr);
41645960306SStefano Zampini     ierr = VecAssemblyEnd(bt);CHKERRQ(ierr);
41745960306SStefano Zampini     ierr = PetscFree(vals);CHKERRQ(ierr);
41845960306SStefano Zampini 
41945960306SStefano Zampini     ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr);
42045960306SStefano Zampini     ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr);
42145960306SStefano Zampini     ierr = VecGetArray(xt,&vals);CHKERRQ(ierr);
42245960306SStefano Zampini     for (i = 0; i < nl; i++) {
42345960306SStefano Zampini       PetscInt g = i + st;
42445960306SStefano Zampini       if (g > mat->rmap->N) continue;
42545960306SStefano Zampini       if (PetscAbsScalar(vals[i]) == 0.0) continue;
42645960306SStefano Zampini       ierr = VecSetValue(b,g,-diag*vals[i],INSERT_VALUES);CHKERRQ(ierr);
42745960306SStefano Zampini     }
42845960306SStefano Zampini     ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr);
42945960306SStefano Zampini     ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
43045960306SStefano Zampini     ierr = VecAssemblyEnd(b);CHKERRQ(ierr);                            /* b  = [b1 - A12*x2, x2 * diag] */
43145960306SStefano Zampini     ierr = VecDestroy(&xt);CHKERRQ(ierr);
43245960306SStefano Zampini     ierr = VecDestroy(&bt);CHKERRQ(ierr);
43345960306SStefano Zampini     ierr = PetscFree2(grows,gcols);CHKERRQ(ierr);
43445960306SStefano Zampini   }
43545960306SStefano Zampini   ierr = PetscLayoutMapLocal(mat->rmap,n,rowscols,&nr,&lrows,NULL);CHKERRQ(ierr);
43645960306SStefano Zampini   ierr = MatHasCongruentLayouts(mat,&congruent);CHKERRQ(ierr);
43745960306SStefano Zampini   if (congruent) {
43845960306SStefano Zampini     nc    = nr;
43945960306SStefano Zampini     lcols = lrows;
44045960306SStefano Zampini   } else { /* MatZeroRowsColumns implicitly assumes the rowscols indices are for a square matrix, here we handle a more general case */
44145960306SStefano Zampini     PetscInt i,nt,*t;
44245960306SStefano Zampini 
44345960306SStefano Zampini     ierr = PetscMalloc1(n,&t);CHKERRQ(ierr);
44445960306SStefano Zampini     for (i = 0, nt = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) t[nt++] = rowscols[i];
44545960306SStefano Zampini     ierr = PetscLayoutMapLocal(mat->cmap,nt,t,&nc,&lcols,NULL);CHKERRQ(ierr);
44645960306SStefano Zampini     ierr = PetscFree(t);CHKERRQ(ierr);
44745960306SStefano Zampini   }
44845960306SStefano Zampini   ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,nc,lcols,diag,PETSC_TRUE);CHKERRQ(ierr);
44945960306SStefano Zampini   if (!congruent) {
45045960306SStefano Zampini     ierr = PetscFree(lcols);CHKERRQ(ierr);
45145960306SStefano Zampini   }
45245960306SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
453*ef5c7bd2SStefano Zampini   if (shell->axpy) {
454*ef5c7bd2SStefano Zampini     ierr = MatZeroRowsColumns(shell->axpy,n,rowscols,0.0,NULL,NULL);CHKERRQ(ierr);
455*ef5c7bd2SStefano Zampini   }
45645960306SStefano Zampini   PetscFunctionReturn(0);
45745960306SStefano Zampini }
45845960306SStefano Zampini 
459dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
460e51e0e81SBarry Smith {
461dfbe8321SBarry Smith   PetscErrorCode ierr;
462bf0cc555SLisandro Dalcin   Mat_Shell      *shell = (Mat_Shell*)mat->data;
463ed3cc1f0SBarry Smith 
4643a40ed3dSBarry Smith   PetscFunctionBegin;
46528f357e3SAlex Fikl   if (shell->ops->destroy) {
46628f357e3SAlex Fikl     ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr);
467bf0cc555SLisandro Dalcin   }
468*ef5c7bd2SStefano Zampini   ierr = PetscMemzero(shell->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
4690c0fd78eSBarry Smith   ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
4700c0fd78eSBarry Smith   ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
4710c0fd78eSBarry Smith   ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
4728fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
4738fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
4745edf6516SJed Brown   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
4755edf6516SJed Brown   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
4769f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
47745960306SStefano Zampini   ierr = VecDestroy(&shell->zvals_w);CHKERRQ(ierr);
47845960306SStefano Zampini   ierr = VecDestroy(&shell->zvals);CHKERRQ(ierr);
47945960306SStefano Zampini   ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr);
48045960306SStefano Zampini   ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr);
48145960306SStefano Zampini   ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr);
48245960306SStefano Zampini   ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr);
483bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
484789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellGetContext_C",NULL);CHKERRQ(ierr);
485789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetContext_C",NULL);CHKERRQ(ierr);
486db77db73SJeremy L Thompson   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetVecType_C",NULL);CHKERRQ(ierr);
487789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetManageScalingShifts_C",NULL);CHKERRQ(ierr);
488789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetOperation_C",NULL);CHKERRQ(ierr);
489789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellGetOperation_C",NULL);CHKERRQ(ierr);
4903a40ed3dSBarry Smith   PetscFunctionReturn(0);
491e51e0e81SBarry Smith }
492e51e0e81SBarry Smith 
4937fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
4947fabda3fSAlex Fikl {
49528f357e3SAlex Fikl   Mat_Shell       *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
4967fabda3fSAlex Fikl   PetscErrorCode  ierr;
497cb8c8a77SBarry Smith   PetscBool       matflg;
4987fabda3fSAlex Fikl 
4997fabda3fSAlex Fikl   PetscFunctionBegin;
50028f357e3SAlex Fikl   ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr);
501cb8c8a77SBarry Smith   if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL");
5027fabda3fSAlex Fikl 
503cb8c8a77SBarry Smith   ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
504cb8c8a77SBarry Smith   ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
5057fabda3fSAlex Fikl 
506cb8c8a77SBarry Smith   if (shellA->ops->copy) {
50728f357e3SAlex Fikl     ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr);
508cb8c8a77SBarry Smith   }
5097fabda3fSAlex Fikl   shellB->vscale = shellA->vscale;
5107fabda3fSAlex Fikl   shellB->vshift = shellA->vshift;
5110c0fd78eSBarry Smith   if (shellA->dshift) {
5120c0fd78eSBarry Smith     if (!shellB->dshift) {
5130c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr);
5147fabda3fSAlex Fikl     }
5150c0fd78eSBarry Smith     ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr);
5167fabda3fSAlex Fikl   } else {
5179f137db4SBarry Smith     ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr);
5187fabda3fSAlex Fikl   }
5190c0fd78eSBarry Smith   if (shellA->left) {
5200c0fd78eSBarry Smith     if (!shellB->left) {
5210c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr);
5227fabda3fSAlex Fikl     }
5230c0fd78eSBarry Smith     ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr);
5247fabda3fSAlex Fikl   } else {
5259f137db4SBarry Smith     ierr = VecDestroy(&shellB->left);CHKERRQ(ierr);
5267fabda3fSAlex Fikl   }
5270c0fd78eSBarry Smith   if (shellA->right) {
5280c0fd78eSBarry Smith     if (!shellB->right) {
5290c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr);
5307fabda3fSAlex Fikl     }
5310c0fd78eSBarry Smith     ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr);
5327fabda3fSAlex Fikl   } else {
5339f137db4SBarry Smith     ierr = VecDestroy(&shellB->right);CHKERRQ(ierr);
5347fabda3fSAlex Fikl   }
53540e381c3SBarry Smith   ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr);
536*ef5c7bd2SStefano Zampini   shellB->axpy_vscale = 0.0;
537*ef5c7bd2SStefano Zampini   shellB->axpy_state  = 0;
53840e381c3SBarry Smith   if (shellA->axpy) {
53940e381c3SBarry Smith     ierr                = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr);
54040e381c3SBarry Smith     shellB->axpy        = shellA->axpy;
54140e381c3SBarry Smith     shellB->axpy_vscale = shellA->axpy_vscale;
542*ef5c7bd2SStefano Zampini     shellB->axpy_state  = shellA->axpy_state;
54340e381c3SBarry Smith   }
54445960306SStefano Zampini   if (shellA->zrows) {
54545960306SStefano Zampini     ierr = ISDuplicate(shellA->zrows,&shellB->zrows);CHKERRQ(ierr);
54645960306SStefano Zampini     if (shellA->zcols) {
54745960306SStefano Zampini       ierr = ISDuplicate(shellA->zcols,&shellB->zcols);CHKERRQ(ierr);
54845960306SStefano Zampini     }
54945960306SStefano Zampini     ierr = VecDuplicate(shellA->zvals,&shellB->zvals);CHKERRQ(ierr);
55045960306SStefano Zampini     ierr = VecCopy(shellA->zvals,shellB->zvals);CHKERRQ(ierr);
55145960306SStefano Zampini     ierr = VecDuplicate(shellA->zvals_w,&shellB->zvals_w);CHKERRQ(ierr);
55245960306SStefano Zampini     ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_r);CHKERRQ(ierr);
55345960306SStefano Zampini     ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_c);CHKERRQ(ierr);
55445960306SStefano Zampini     shellB->zvals_sct_r = shellA->zvals_sct_r;
55545960306SStefano Zampini     shellB->zvals_sct_c = shellA->zvals_sct_c;
55645960306SStefano Zampini   }
5577fabda3fSAlex Fikl   PetscFunctionReturn(0);
5587fabda3fSAlex Fikl }
5597fabda3fSAlex Fikl 
560cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M)
561cb8c8a77SBarry Smith {
562cb8c8a77SBarry Smith   PetscErrorCode ierr;
563cb8c8a77SBarry Smith   void           *ctx;
564cb8c8a77SBarry Smith 
565cb8c8a77SBarry Smith   PetscFunctionBegin;
566cb8c8a77SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
567cb8c8a77SBarry Smith   ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr);
568a4b1380bSStefano Zampini   if (op != MAT_DO_NOT_COPY_VALUES) {
569cb8c8a77SBarry Smith     ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
570a4b1380bSStefano Zampini   }
571cb8c8a77SBarry Smith   PetscFunctionReturn(0);
572cb8c8a77SBarry Smith }
573cb8c8a77SBarry Smith 
574dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
575ef66eb69SBarry Smith {
576ef66eb69SBarry Smith   Mat_Shell        *shell = (Mat_Shell*)A->data;
577dfbe8321SBarry Smith   PetscErrorCode   ierr;
57825578ef6SJed Brown   Vec              xx;
579e3f487b0SBarry Smith   PetscObjectState instate,outstate;
580ef66eb69SBarry Smith 
581ef66eb69SBarry Smith   PetscFunctionBegin;
582976e8c5aSLisandro Dalcin   if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL");
58345960306SStefano Zampini   ierr = MatShellPreZeroRight(A,x,&xx);CHKERRQ(ierr);
58445960306SStefano Zampini   ierr = MatShellPreScaleRight(A,xx,&xx);CHKERRQ(ierr);
58545960306SStefano Zampini   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
58628f357e3SAlex Fikl   ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr);
587e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
588e3f487b0SBarry Smith   if (instate == outstate) {
589e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
590e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
591e3f487b0SBarry Smith   }
5928fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
5938fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
59445960306SStefano Zampini   ierr = MatShellPostZeroLeft(A,y);CHKERRQ(ierr);
5959f137db4SBarry Smith 
5969f137db4SBarry Smith   if (shell->axpy) {
597*ef5c7bd2SStefano Zampini     Mat              X;
598*ef5c7bd2SStefano Zampini     PetscObjectState axpy_state;
599*ef5c7bd2SStefano Zampini 
600*ef5c7bd2SStefano Zampini     ierr = MatShellGetContext(shell->axpy,(void *)&X);CHKERRQ(ierr);
601*ef5c7bd2SStefano Zampini     ierr = PetscObjectStateGet((PetscObject)X,&axpy_state);CHKERRQ(ierr);
602*ef5c7bd2SStefano Zampini     if (shell->axpy_state != axpy_state) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state %D != %D",axpy_state,shell->axpy_state);
6039f137db4SBarry Smith     if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);}
6049f137db4SBarry Smith     ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr);
6059f137db4SBarry Smith     ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
6069f137db4SBarry Smith   }
607ef66eb69SBarry Smith   PetscFunctionReturn(0);
608ef66eb69SBarry Smith }
609ef66eb69SBarry Smith 
6105edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
6115edf6516SJed Brown {
6125edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
6135edf6516SJed Brown   PetscErrorCode ierr;
6145edf6516SJed Brown 
6155edf6516SJed Brown   PetscFunctionBegin;
6165edf6516SJed Brown   if (y == z) {
6175edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
6185edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
619b8430d49SBarry Smith     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
6205edf6516SJed Brown   } else {
6215edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
6225edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
6235edf6516SJed Brown   }
6245edf6516SJed Brown   PetscFunctionReturn(0);
6255edf6516SJed Brown }
6265edf6516SJed Brown 
62718be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
62818be62a5SSatish Balay {
62918be62a5SSatish Balay   Mat_Shell        *shell = (Mat_Shell*)A->data;
63018be62a5SSatish Balay   PetscErrorCode   ierr;
63125578ef6SJed Brown   Vec              xx;
632e3f487b0SBarry Smith   PetscObjectState instate,outstate;
63318be62a5SSatish Balay 
63418be62a5SSatish Balay   PetscFunctionBegin;
63545960306SStefano Zampini   ierr = MatShellPreZeroLeft(A,x,&xx);CHKERRQ(ierr);
63645960306SStefano Zampini   ierr = MatShellPreScaleLeft(A,xx,&xx);CHKERRQ(ierr);
637e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
6380c0fd78eSBarry Smith   if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL");
63928f357e3SAlex Fikl   ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr);
640e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
641e3f487b0SBarry Smith   if (instate == outstate) {
642e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
643e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
644e3f487b0SBarry Smith   }
6458fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
6468fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
64745960306SStefano Zampini   ierr = MatShellPostZeroRight(A,y);CHKERRQ(ierr);
648350ec83bSStefano Zampini 
649350ec83bSStefano Zampini   if (shell->axpy) {
650*ef5c7bd2SStefano Zampini     Mat              X;
651*ef5c7bd2SStefano Zampini     PetscObjectState axpy_state;
652*ef5c7bd2SStefano Zampini 
653*ef5c7bd2SStefano Zampini     ierr = MatShellGetContext(shell->axpy,(void *)&X);CHKERRQ(ierr);
654*ef5c7bd2SStefano Zampini     ierr = PetscObjectStateGet((PetscObject)X,&axpy_state);CHKERRQ(ierr);
655*ef5c7bd2SStefano Zampini     if (shell->axpy_state != axpy_state) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state %D != %D",axpy_state,shell->axpy_state);
656350ec83bSStefano Zampini     if (!shell->right_work) {ierr = MatCreateVecs(A,NULL,&shell->right_work);CHKERRQ(ierr);}
657350ec83bSStefano Zampini     ierr = MatMultTranspose(shell->axpy,x,shell->right_work);CHKERRQ(ierr);
658350ec83bSStefano Zampini     ierr = VecAXPY(y,shell->axpy_vscale,shell->right_work);CHKERRQ(ierr);
659350ec83bSStefano Zampini   }
66018be62a5SSatish Balay   PetscFunctionReturn(0);
66118be62a5SSatish Balay }
66218be62a5SSatish Balay 
6635edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
6645edf6516SJed Brown {
6655edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
6665edf6516SJed Brown   PetscErrorCode ierr;
6675edf6516SJed Brown 
6685edf6516SJed Brown   PetscFunctionBegin;
6695edf6516SJed Brown   if (y == z) {
6705edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
6715edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
672dac989eaS“Marek     ierr = VecAXPY(z,1.0,shell->left_add_work);CHKERRQ(ierr);
6735edf6516SJed Brown   } else {
6745edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
6755edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
6765edf6516SJed Brown   }
6775edf6516SJed Brown   PetscFunctionReturn(0);
6785edf6516SJed Brown }
6795edf6516SJed Brown 
6800c0fd78eSBarry Smith /*
6810c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
6820c0fd78eSBarry Smith */
68318be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
68418be62a5SSatish Balay {
68518be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
68618be62a5SSatish Balay   PetscErrorCode ierr;
68718be62a5SSatish Balay 
68818be62a5SSatish Balay   PetscFunctionBegin;
6890c0fd78eSBarry Smith   if (shell->ops->getdiagonal) {
69028f357e3SAlex Fikl     ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr);
691305211d5SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Must provide shell matrix with routine to return diagonal using\nMatShellSetOperation(S,MATOP_GET_DIAGONAL,...)");
69218be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
6938fe8eb27SJed Brown   if (shell->dshift) {
6940c0fd78eSBarry Smith     ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr);
69518be62a5SSatish Balay   }
6960c0fd78eSBarry Smith   ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
6978fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
6988fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
69945960306SStefano Zampini   if (shell->zrows) {
70045960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
70145960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
70245960306SStefano Zampini   }
70381c02519SBarry Smith   if (shell->axpy) {
704*ef5c7bd2SStefano Zampini     Mat              X;
705*ef5c7bd2SStefano Zampini     PetscObjectState axpy_state;
706*ef5c7bd2SStefano Zampini 
707*ef5c7bd2SStefano Zampini     ierr = MatShellGetContext(shell->axpy,(void *)&X);CHKERRQ(ierr);
708*ef5c7bd2SStefano Zampini     ierr = PetscObjectStateGet((PetscObject)X,&axpy_state);CHKERRQ(ierr);
709*ef5c7bd2SStefano Zampini     if (shell->axpy_state != axpy_state) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state %D != %D",axpy_state,shell->axpy_state);
71081c02519SBarry Smith     if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);}
71181c02519SBarry Smith     ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr);
71281c02519SBarry Smith     ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
71381c02519SBarry Smith   }
71418be62a5SSatish Balay   PetscFunctionReturn(0);
71518be62a5SSatish Balay }
71618be62a5SSatish Balay 
717f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
718ef66eb69SBarry Smith {
719ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
7208fe8eb27SJed Brown   PetscErrorCode ierr;
721789d8953SBarry Smith   PetscBool      flg;
722b24ad042SBarry Smith 
723ef66eb69SBarry Smith   PetscFunctionBegin;
724789d8953SBarry Smith   ierr = MatHasCongruentLayouts(Y,&flg);CHKERRQ(ierr);
725789d8953SBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot shift shell matrix if it is not congruent");
7260c0fd78eSBarry Smith   if (shell->left || shell->right) {
7278fe8eb27SJed Brown     if (!shell->dshift) {
7280c0fd78eSBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
7290c0fd78eSBarry Smith       ierr = VecSet(shell->dshift,a);CHKERRQ(ierr);
7300c0fd78eSBarry Smith     } else {
7310c0fd78eSBarry Smith       if (shell->left)  {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
7320c0fd78eSBarry Smith       if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
7339f137db4SBarry Smith       ierr = VecShift(shell->dshift,a);CHKERRQ(ierr);
7340c0fd78eSBarry Smith     }
7358fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
7368fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
7378fe8eb27SJed Brown   } else shell->vshift += a;
73845960306SStefano Zampini   if (shell->zrows) {
73945960306SStefano Zampini     ierr = VecShift(shell->zvals,a);CHKERRQ(ierr);
74045960306SStefano Zampini   }
741ef66eb69SBarry Smith   PetscFunctionReturn(0);
742ef66eb69SBarry Smith }
743ef66eb69SBarry Smith 
744b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s)
7456464bf51SAlex Fikl {
7466464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
7476464bf51SAlex Fikl   PetscErrorCode ierr;
7486464bf51SAlex Fikl 
7496464bf51SAlex Fikl   PetscFunctionBegin;
7500c0fd78eSBarry Smith   if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);}
7510c0fd78eSBarry Smith   if (shell->left || shell->right) {
75292fabd1fSBarry Smith     if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);}
7539f137db4SBarry Smith     if (shell->left && shell->right)  {
7549f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
7559f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
7569f137db4SBarry Smith     } else if (shell->left) {
7579f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
7589f137db4SBarry Smith     } else {
7599f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr);
7609f137db4SBarry Smith     }
761b253ae0bS“Marek     ierr = VecAXPY(shell->dshift,s,shell->right_work);CHKERRQ(ierr);
7620c0fd78eSBarry Smith   } else {
763b253ae0bS“Marek     ierr = VecAXPY(shell->dshift,s,D);CHKERRQ(ierr);
764b253ae0bS“Marek   }
765b253ae0bS“Marek   PetscFunctionReturn(0);
766b253ae0bS“Marek }
767b253ae0bS“Marek 
768b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
769b253ae0bS“Marek {
77045960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
771b253ae0bS“Marek   Vec            d;
772b253ae0bS“Marek   PetscErrorCode ierr;
773789d8953SBarry Smith   PetscBool      flg;
774b253ae0bS“Marek 
775b253ae0bS“Marek   PetscFunctionBegin;
776789d8953SBarry Smith   ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr);
777789d8953SBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot diagonal set or shift shell matrix if it is not congruent");
778b253ae0bS“Marek   if (ins == INSERT_VALUES) {
779b253ae0bS“Marek     if (!A->ops->getdiagonal) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first");
780b253ae0bS“Marek     ierr = VecDuplicate(D,&d);CHKERRQ(ierr);
781b253ae0bS“Marek     ierr = MatGetDiagonal(A,d);CHKERRQ(ierr);
782b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,d,-1.);CHKERRQ(ierr);
783b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr);
784b253ae0bS“Marek     ierr = VecDestroy(&d);CHKERRQ(ierr);
78545960306SStefano Zampini     if (shell->zrows) {
78645960306SStefano Zampini       ierr = VecCopy(D,shell->zvals);CHKERRQ(ierr);
78745960306SStefano Zampini     }
788b253ae0bS“Marek   } else {
789b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr);
79045960306SStefano Zampini     if (shell->zrows) {
79145960306SStefano Zampini       ierr = VecAXPY(shell->zvals,1.0,D);CHKERRQ(ierr);
79245960306SStefano Zampini     }
7936464bf51SAlex Fikl   }
7946464bf51SAlex Fikl   PetscFunctionReturn(0);
7956464bf51SAlex Fikl }
7966464bf51SAlex Fikl 
797f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
798ef66eb69SBarry Smith {
799ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
8008fe8eb27SJed Brown   PetscErrorCode ierr;
801b24ad042SBarry Smith 
802ef66eb69SBarry Smith   PetscFunctionBegin;
803f4df32b1SMatthew Knepley   shell->vscale *= a;
8040c0fd78eSBarry Smith   shell->vshift *= a;
8052205254eSKarl Rupp   if (shell->dshift) {
8062205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
8070c0fd78eSBarry Smith   }
80881c02519SBarry Smith   shell->axpy_vscale *= a;
80945960306SStefano Zampini   if (shell->zrows) {
81045960306SStefano Zampini     ierr = VecScale(shell->zvals,a);CHKERRQ(ierr);
81145960306SStefano Zampini   }
8128fe8eb27SJed Brown   PetscFunctionReturn(0);
81318be62a5SSatish Balay }
8148fe8eb27SJed Brown 
8158fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
8168fe8eb27SJed Brown {
8178fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
8188fe8eb27SJed Brown   PetscErrorCode ierr;
8198fe8eb27SJed Brown 
8208fe8eb27SJed Brown   PetscFunctionBegin;
8218fe8eb27SJed Brown   if (left) {
8220c0fd78eSBarry Smith     if (!shell->left) {
8230c0fd78eSBarry Smith       ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr);
8248fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
8250c0fd78eSBarry Smith     } else {
8260c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
82718be62a5SSatish Balay     }
82845960306SStefano Zampini     if (shell->zrows) {
82945960306SStefano Zampini       ierr = VecPointwiseMult(shell->zvals,shell->zvals,left);CHKERRQ(ierr);
83045960306SStefano Zampini     }
831ef66eb69SBarry Smith   }
8328fe8eb27SJed Brown   if (right) {
8330c0fd78eSBarry Smith     if (!shell->right) {
8340c0fd78eSBarry Smith       ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr);
8358fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
8360c0fd78eSBarry Smith     } else {
8370c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
8388fe8eb27SJed Brown     }
83945960306SStefano Zampini     if (shell->zrows) {
84045960306SStefano Zampini       if (!shell->left_work) {
84145960306SStefano Zampini         ierr = MatCreateVecs(Y,NULL,&shell->left_work);CHKERRQ(ierr);
84245960306SStefano Zampini       }
84345960306SStefano Zampini       ierr = VecSet(shell->zvals_w,1.0);CHKERRQ(ierr);
84445960306SStefano Zampini       ierr = VecScatterBegin(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
84545960306SStefano Zampini       ierr = VecScatterEnd(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
84645960306SStefano Zampini       ierr = VecPointwiseMult(shell->zvals,shell->zvals,shell->zvals_w);CHKERRQ(ierr);
84745960306SStefano Zampini     }
8488fe8eb27SJed Brown   }
849*ef5c7bd2SStefano Zampini   if (shell->axpy) {
850*ef5c7bd2SStefano Zampini     ierr = MatDiagonalScale(shell->axpy,left,right);CHKERRQ(ierr);
851*ef5c7bd2SStefano Zampini   }
852ef66eb69SBarry Smith   PetscFunctionReturn(0);
853ef66eb69SBarry Smith }
854ef66eb69SBarry Smith 
855dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
856ef66eb69SBarry Smith {
857ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
8580c0fd78eSBarry Smith   PetscErrorCode ierr;
859ef66eb69SBarry Smith 
860ef66eb69SBarry Smith   PetscFunctionBegin;
8618fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
862ef66eb69SBarry Smith     shell->vshift = 0.0;
863ef66eb69SBarry Smith     shell->vscale = 1.0;
864*ef5c7bd2SStefano Zampini     shell->axpy_vscale = 0.0;
865*ef5c7bd2SStefano Zampini     shell->axpy_state  = 0;
8660c0fd78eSBarry Smith     ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
8670c0fd78eSBarry Smith     ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
8680c0fd78eSBarry Smith     ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
86981c02519SBarry Smith     ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
87045960306SStefano Zampini     ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr);
87145960306SStefano Zampini     ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr);
87245960306SStefano Zampini     ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr);
87345960306SStefano Zampini     ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr);
874ef66eb69SBarry Smith   }
875ef66eb69SBarry Smith   PetscFunctionReturn(0);
876ef66eb69SBarry Smith }
877ef66eb69SBarry Smith 
8783b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d)
8793b49f96aSBarry Smith {
8803b49f96aSBarry Smith   PetscFunctionBegin;
8813b49f96aSBarry Smith   *missing = PETSC_FALSE;
8823b49f96aSBarry Smith   PetscFunctionReturn(0);
8833b49f96aSBarry Smith }
8843b49f96aSBarry Smith 
8859f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
8869f137db4SBarry Smith {
8879f137db4SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
8889f137db4SBarry Smith   PetscErrorCode ierr;
8899f137db4SBarry Smith 
8909f137db4SBarry Smith   PetscFunctionBegin;
891*ef5c7bd2SStefano Zampini   if (X == Y) {
892*ef5c7bd2SStefano Zampini     ierr = MatScale(Y,1.0 + a);CHKERRQ(ierr);
893*ef5c7bd2SStefano Zampini     PetscFunctionReturn(0);
894*ef5c7bd2SStefano Zampini   }
895*ef5c7bd2SStefano Zampini   if (!shell->axpy) {
896*ef5c7bd2SStefano Zampini     ierr = MatConvertFrom_Shell(X,MATSHELL,MAT_INITIAL_MATRIX,&shell->axpy);CHKERRQ(ierr);
8979f137db4SBarry Smith     shell->axpy_vscale = a;
898*ef5c7bd2SStefano Zampini     ierr = PetscObjectStateGet((PetscObject)X,&shell->axpy_state);CHKERRQ(ierr);
899*ef5c7bd2SStefano Zampini   } else {
900*ef5c7bd2SStefano Zampini     ierr = MatAXPY(shell->axpy,a/shell->axpy_vscale,X,str);CHKERRQ(ierr);
901*ef5c7bd2SStefano Zampini   }
9029f137db4SBarry Smith   PetscFunctionReturn(0);
9039f137db4SBarry Smith }
9049f137db4SBarry Smith 
90509dc0095SBarry Smith static struct _MatOps MatOps_Values = {0,
90620563c6bSBarry Smith                                        0,
90720563c6bSBarry Smith                                        0,
9089f137db4SBarry Smith                                        0,
9090c0fd78eSBarry Smith                                 /* 4*/ MatMultAdd_Shell,
9109f137db4SBarry Smith                                        0,
9110c0fd78eSBarry Smith                                        MatMultTransposeAdd_Shell,
912b951964fSBarry Smith                                        0,
913b951964fSBarry Smith                                        0,
914b951964fSBarry Smith                                        0,
91597304618SKris Buschelman                                 /*10*/ 0,
916b951964fSBarry Smith                                        0,
917b951964fSBarry Smith                                        0,
918b951964fSBarry Smith                                        0,
919b951964fSBarry Smith                                        0,
92097304618SKris Buschelman                                 /*15*/ 0,
921b951964fSBarry Smith                                        0,
92200849d43SBarry Smith                                        0,
9238fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
924b951964fSBarry Smith                                        0,
92597304618SKris Buschelman                                 /*20*/ 0,
926ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
927b951964fSBarry Smith                                        0,
928b951964fSBarry Smith                                        0,
92945960306SStefano Zampini                                 /*24*/ MatZeroRows_Shell,
930b951964fSBarry Smith                                        0,
931b951964fSBarry Smith                                        0,
932b951964fSBarry Smith                                        0,
933b951964fSBarry Smith                                        0,
934d519adbfSMatthew Knepley                                 /*29*/ 0,
935b951964fSBarry Smith                                        0,
936273d9f13SBarry Smith                                        0,
937b951964fSBarry Smith                                        0,
938b951964fSBarry Smith                                        0,
939cb8c8a77SBarry Smith                                 /*34*/ MatDuplicate_Shell,
940b951964fSBarry Smith                                        0,
941b951964fSBarry Smith                                        0,
94209dc0095SBarry Smith                                        0,
94309dc0095SBarry Smith                                        0,
9449f137db4SBarry Smith                                 /*39*/ MatAXPY_Shell,
94509dc0095SBarry Smith                                        0,
94609dc0095SBarry Smith                                        0,
94709dc0095SBarry Smith                                        0,
948cb8c8a77SBarry Smith                                        MatCopy_Shell,
949d519adbfSMatthew Knepley                                 /*44*/ 0,
950ef66eb69SBarry Smith                                        MatScale_Shell,
951ef66eb69SBarry Smith                                        MatShift_Shell,
9526464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
95345960306SStefano Zampini                                        MatZeroRowsColumns_Shell,
954f73d5cc4SBarry Smith                                 /*49*/ 0,
95509dc0095SBarry Smith                                        0,
95609dc0095SBarry Smith                                        0,
95709dc0095SBarry Smith                                        0,
95809dc0095SBarry Smith                                        0,
959d519adbfSMatthew Knepley                                 /*54*/ 0,
96009dc0095SBarry Smith                                        0,
96109dc0095SBarry Smith                                        0,
96209dc0095SBarry Smith                                        0,
96309dc0095SBarry Smith                                        0,
964d519adbfSMatthew Knepley                                 /*59*/ 0,
965b9b97703SBarry Smith                                        MatDestroy_Shell,
96609dc0095SBarry Smith                                        0,
967251fad3fSStefano Zampini                                        MatConvertFrom_Shell,
968273d9f13SBarry Smith                                        0,
969d519adbfSMatthew Knepley                                 /*64*/ 0,
970273d9f13SBarry Smith                                        0,
971273d9f13SBarry Smith                                        0,
972273d9f13SBarry Smith                                        0,
973273d9f13SBarry Smith                                        0,
974d519adbfSMatthew Knepley                                 /*69*/ 0,
975273d9f13SBarry Smith                                        0,
976c87e5d42SMatthew Knepley                                        MatConvert_Shell,
977273d9f13SBarry Smith                                        0,
97897304618SKris Buschelman                                        0,
979d519adbfSMatthew Knepley                                 /*74*/ 0,
98097304618SKris Buschelman                                        0,
98197304618SKris Buschelman                                        0,
98297304618SKris Buschelman                                        0,
98397304618SKris Buschelman                                        0,
984d519adbfSMatthew Knepley                                 /*79*/ 0,
98597304618SKris Buschelman                                        0,
98697304618SKris Buschelman                                        0,
98797304618SKris Buschelman                                        0,
988865e5f61SKris Buschelman                                        0,
989d519adbfSMatthew Knepley                                 /*84*/ 0,
990865e5f61SKris Buschelman                                        0,
991865e5f61SKris Buschelman                                        0,
992865e5f61SKris Buschelman                                        0,
993865e5f61SKris Buschelman                                        0,
994d519adbfSMatthew Knepley                                 /*89*/ 0,
995865e5f61SKris Buschelman                                        0,
996865e5f61SKris Buschelman                                        0,
997865e5f61SKris Buschelman                                        0,
998865e5f61SKris Buschelman                                        0,
999d519adbfSMatthew Knepley                                 /*94*/ 0,
1000865e5f61SKris Buschelman                                        0,
1001865e5f61SKris Buschelman                                        0,
10023964eb88SJed Brown                                        0,
10033964eb88SJed Brown                                        0,
10043964eb88SJed Brown                                 /*99*/ 0,
10053964eb88SJed Brown                                        0,
10063964eb88SJed Brown                                        0,
10073964eb88SJed Brown                                        0,
10083964eb88SJed Brown                                        0,
10093964eb88SJed Brown                                /*104*/ 0,
10103964eb88SJed Brown                                        0,
10113964eb88SJed Brown                                        0,
10123964eb88SJed Brown                                        0,
10133964eb88SJed Brown                                        0,
10143964eb88SJed Brown                                /*109*/ 0,
10153964eb88SJed Brown                                        0,
10163964eb88SJed Brown                                        0,
10173964eb88SJed Brown                                        0,
10183b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
10193964eb88SJed Brown                                /*114*/ 0,
10203964eb88SJed Brown                                        0,
10213964eb88SJed Brown                                        0,
10223964eb88SJed Brown                                        0,
10233964eb88SJed Brown                                        0,
10243964eb88SJed Brown                                /*119*/ 0,
10253964eb88SJed Brown                                        0,
10263964eb88SJed Brown                                        0,
10273964eb88SJed Brown                                        0,
10283964eb88SJed Brown                                        0,
10293964eb88SJed Brown                                /*124*/ 0,
10303964eb88SJed Brown                                        0,
10313964eb88SJed Brown                                        0,
10323964eb88SJed Brown                                        0,
10333964eb88SJed Brown                                        0,
10343964eb88SJed Brown                                /*129*/ 0,
10353964eb88SJed Brown                                        0,
10363964eb88SJed Brown                                        0,
10373964eb88SJed Brown                                        0,
10383964eb88SJed Brown                                        0,
10393964eb88SJed Brown                                /*134*/ 0,
10403964eb88SJed Brown                                        0,
10413964eb88SJed Brown                                        0,
10423964eb88SJed Brown                                        0,
10433964eb88SJed Brown                                        0,
10443964eb88SJed Brown                                /*139*/ 0,
1045f9426fe0SMark Adams                                        0,
10463964eb88SJed Brown                                        0
10473964eb88SJed Brown };
1048273d9f13SBarry Smith 
1049789d8953SBarry Smith PetscErrorCode  MatShellSetContext_Shell(Mat mat,void *ctx)
1050789d8953SBarry Smith {
1051789d8953SBarry Smith   Mat_Shell *shell = (Mat_Shell*)mat->data;
1052789d8953SBarry Smith 
1053789d8953SBarry Smith   PetscFunctionBegin;
1054789d8953SBarry Smith   shell->ctx = ctx;
1055789d8953SBarry Smith   PetscFunctionReturn(0);
1056789d8953SBarry Smith }
1057789d8953SBarry Smith 
1058db77db73SJeremy L Thompson static PetscErrorCode MatShellSetVecType_Shell(Mat mat,VecType vtype)
1059db77db73SJeremy L Thompson {
1060db77db73SJeremy L Thompson   PetscErrorCode ierr;
1061db77db73SJeremy L Thompson 
1062db77db73SJeremy L Thompson   PetscFunctionBegin;
1063db77db73SJeremy L Thompson   ierr = PetscFree(mat->defaultvectype);CHKERRQ(ierr);
1064db77db73SJeremy L Thompson   ierr = PetscStrallocpy(vtype,(char**)&mat->defaultvectype);CHKERRQ(ierr);
1065db77db73SJeremy L Thompson   PetscFunctionReturn(0);
1066db77db73SJeremy L Thompson }
1067db77db73SJeremy L Thompson 
1068789d8953SBarry Smith PetscErrorCode MatShellSetManageScalingShifts_Shell(Mat A)
1069789d8953SBarry Smith {
1070789d8953SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)A->data;
1071789d8953SBarry Smith 
1072789d8953SBarry Smith   PetscFunctionBegin;
1073789d8953SBarry Smith   shell->managescalingshifts = PETSC_FALSE;
1074789d8953SBarry Smith   A->ops->diagonalset   = NULL;
1075789d8953SBarry Smith   A->ops->diagonalscale = NULL;
1076789d8953SBarry Smith   A->ops->scale         = NULL;
1077789d8953SBarry Smith   A->ops->shift         = NULL;
1078789d8953SBarry Smith   A->ops->axpy          = NULL;
1079789d8953SBarry Smith   PetscFunctionReturn(0);
1080789d8953SBarry Smith }
1081789d8953SBarry Smith 
1082789d8953SBarry Smith PetscErrorCode MatShellSetOperation_Shell(Mat mat,MatOperation op,void (*f)(void))
1083789d8953SBarry Smith {
1084feb237baSPierre Jolivet   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1085789d8953SBarry Smith 
1086789d8953SBarry Smith   PetscFunctionBegin;
1087789d8953SBarry Smith   switch (op) {
1088789d8953SBarry Smith   case MATOP_DESTROY:
1089789d8953SBarry Smith     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
1090789d8953SBarry Smith     break;
1091789d8953SBarry Smith   case MATOP_VIEW:
1092789d8953SBarry Smith     if (!mat->ops->viewnative) {
1093789d8953SBarry Smith       mat->ops->viewnative = mat->ops->view;
1094789d8953SBarry Smith     }
1095789d8953SBarry Smith     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
1096789d8953SBarry Smith     break;
1097789d8953SBarry Smith   case MATOP_COPY:
1098789d8953SBarry Smith     shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f;
1099789d8953SBarry Smith     break;
1100789d8953SBarry Smith   case MATOP_DIAGONAL_SET:
1101789d8953SBarry Smith   case MATOP_DIAGONAL_SCALE:
1102789d8953SBarry Smith   case MATOP_SHIFT:
1103789d8953SBarry Smith   case MATOP_SCALE:
1104789d8953SBarry Smith   case MATOP_AXPY:
1105789d8953SBarry Smith   case MATOP_ZERO_ROWS:
1106789d8953SBarry Smith   case MATOP_ZERO_ROWS_COLUMNS:
1107789d8953SBarry Smith     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
1108789d8953SBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
1109789d8953SBarry Smith     break;
1110789d8953SBarry Smith   case MATOP_GET_DIAGONAL:
1111789d8953SBarry Smith     if (shell->managescalingshifts) {
1112789d8953SBarry Smith       shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1113789d8953SBarry Smith       mat->ops->getdiagonal   = MatGetDiagonal_Shell;
1114789d8953SBarry Smith     } else {
1115789d8953SBarry Smith       shell->ops->getdiagonal = NULL;
1116789d8953SBarry Smith       mat->ops->getdiagonal   = (PetscErrorCode (*)(Mat,Vec))f;
1117789d8953SBarry Smith     }
1118789d8953SBarry Smith     break;
1119789d8953SBarry Smith   case MATOP_MULT:
1120789d8953SBarry Smith     if (shell->managescalingshifts) {
1121789d8953SBarry Smith       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1122789d8953SBarry Smith       mat->ops->mult   = MatMult_Shell;
1123789d8953SBarry Smith     } else {
1124789d8953SBarry Smith       shell->ops->mult = NULL;
1125789d8953SBarry Smith       mat->ops->mult   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1126789d8953SBarry Smith     }
1127789d8953SBarry Smith     break;
1128789d8953SBarry Smith   case MATOP_MULT_TRANSPOSE:
1129789d8953SBarry Smith     if (shell->managescalingshifts) {
1130789d8953SBarry Smith       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1131789d8953SBarry Smith       mat->ops->multtranspose   = MatMultTranspose_Shell;
1132789d8953SBarry Smith     } else {
1133789d8953SBarry Smith       shell->ops->multtranspose = NULL;
1134789d8953SBarry Smith       mat->ops->multtranspose   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1135789d8953SBarry Smith     }
1136789d8953SBarry Smith     break;
1137789d8953SBarry Smith   default:
1138789d8953SBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
1139789d8953SBarry Smith     break;
1140789d8953SBarry Smith   }
1141789d8953SBarry Smith   PetscFunctionReturn(0);
1142789d8953SBarry Smith }
1143789d8953SBarry Smith 
1144789d8953SBarry Smith PetscErrorCode MatShellGetOperation_Shell(Mat mat,MatOperation op,void(**f)(void))
1145789d8953SBarry Smith {
1146789d8953SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1147789d8953SBarry Smith 
1148789d8953SBarry Smith   PetscFunctionBegin;
1149789d8953SBarry Smith   switch (op) {
1150789d8953SBarry Smith   case MATOP_DESTROY:
1151789d8953SBarry Smith     *f = (void (*)(void))shell->ops->destroy;
1152789d8953SBarry Smith     break;
1153789d8953SBarry Smith   case MATOP_VIEW:
1154789d8953SBarry Smith     *f = (void (*)(void))mat->ops->view;
1155789d8953SBarry Smith     break;
1156789d8953SBarry Smith   case MATOP_COPY:
1157789d8953SBarry Smith     *f = (void (*)(void))shell->ops->copy;
1158789d8953SBarry Smith     break;
1159789d8953SBarry Smith   case MATOP_DIAGONAL_SET:
1160789d8953SBarry Smith   case MATOP_DIAGONAL_SCALE:
1161789d8953SBarry Smith   case MATOP_SHIFT:
1162789d8953SBarry Smith   case MATOP_SCALE:
1163789d8953SBarry Smith   case MATOP_AXPY:
1164789d8953SBarry Smith   case MATOP_ZERO_ROWS:
1165789d8953SBarry Smith   case MATOP_ZERO_ROWS_COLUMNS:
1166789d8953SBarry Smith     *f = (((void (**)(void))mat->ops)[op]);
1167789d8953SBarry Smith     break;
1168789d8953SBarry Smith   case MATOP_GET_DIAGONAL:
1169789d8953SBarry Smith     if (shell->ops->getdiagonal)
1170789d8953SBarry Smith       *f = (void (*)(void))shell->ops->getdiagonal;
1171789d8953SBarry Smith     else
1172789d8953SBarry Smith       *f = (((void (**)(void))mat->ops)[op]);
1173789d8953SBarry Smith     break;
1174789d8953SBarry Smith   case MATOP_MULT:
1175789d8953SBarry Smith     if (shell->ops->mult)
1176789d8953SBarry Smith       *f = (void (*)(void))shell->ops->mult;
1177789d8953SBarry Smith     else
1178789d8953SBarry Smith       *f = (((void (**)(void))mat->ops)[op]);
1179789d8953SBarry Smith     break;
1180789d8953SBarry Smith   case MATOP_MULT_TRANSPOSE:
1181789d8953SBarry Smith     if (shell->ops->multtranspose)
1182789d8953SBarry Smith       *f = (void (*)(void))shell->ops->multtranspose;
1183789d8953SBarry Smith     else
1184789d8953SBarry Smith       *f = (((void (**)(void))mat->ops)[op]);
1185789d8953SBarry Smith     break;
1186789d8953SBarry Smith   default:
1187789d8953SBarry Smith     *f = (((void (**)(void))mat->ops)[op]);
1188789d8953SBarry Smith   }
1189789d8953SBarry Smith   PetscFunctionReturn(0);
1190789d8953SBarry Smith }
1191789d8953SBarry Smith 
11920bad9183SKris Buschelman /*MC
1193fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
11940bad9183SKris Buschelman 
11950bad9183SKris Buschelman   Level: advanced
11960bad9183SKris Buschelman 
11970c0fd78eSBarry Smith .seealso: MatCreateShell()
11980bad9183SKris Buschelman M*/
11990bad9183SKris Buschelman 
12008cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
1201273d9f13SBarry Smith {
1202273d9f13SBarry Smith   Mat_Shell      *b;
1203dfbe8321SBarry Smith   PetscErrorCode ierr;
1204273d9f13SBarry Smith 
1205273d9f13SBarry Smith   PetscFunctionBegin;
1206273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1207273d9f13SBarry Smith 
1208b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1209273d9f13SBarry Smith   A->data = (void*)b;
1210273d9f13SBarry Smith 
121126283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
121226283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1213273d9f13SBarry Smith 
1214273d9f13SBarry Smith   b->ctx                 = 0;
1215ef66eb69SBarry Smith   b->vshift              = 0.0;
1216ef66eb69SBarry Smith   b->vscale              = 1.0;
12170c0fd78eSBarry Smith   b->managescalingshifts = PETSC_TRUE;
1218273d9f13SBarry Smith   A->assembled           = PETSC_TRUE;
1219210f0121SBarry Smith   A->preallocated        = PETSC_FALSE;
12202205254eSKarl Rupp 
1221789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetContext_C",MatShellGetContext_Shell);CHKERRQ(ierr);
1222789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetContext_C",MatShellSetContext_Shell);CHKERRQ(ierr);
1223db77db73SJeremy L Thompson   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetVecType_C",MatShellSetVecType_Shell);CHKERRQ(ierr);
1224789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetManageScalingShifts_C",MatShellSetManageScalingShifts_Shell);CHKERRQ(ierr);
1225789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetOperation_C",MatShellSetOperation_Shell);CHKERRQ(ierr);
1226789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetOperation_C",MatShellGetOperation_Shell);CHKERRQ(ierr);
122717667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
1228273d9f13SBarry Smith   PetscFunctionReturn(0);
1229273d9f13SBarry Smith }
1230e51e0e81SBarry Smith 
12314b828684SBarry Smith /*@C
1232052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
1233ff756334SLois Curfman McInnes    private data storage format.
1234e51e0e81SBarry Smith 
1235d083f849SBarry Smith   Collective
1236c7fcc2eaSBarry Smith 
1237e51e0e81SBarry Smith    Input Parameters:
1238c7fcc2eaSBarry Smith +  comm - MPI communicator
1239c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
1240c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
1241c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
1242c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
1243c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
1244e51e0e81SBarry Smith 
1245ff756334SLois Curfman McInnes    Output Parameter:
124644cd7ae7SLois Curfman McInnes .  A - the matrix
1247e51e0e81SBarry Smith 
1248ff2fd236SBarry Smith    Level: advanced
1249ff2fd236SBarry Smith 
1250f39d1f56SLois Curfman McInnes   Usage:
12515bd1e576SStefano Zampini $    extern PetscErrorCode mult(Mat,Vec,Vec);
1252f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
1253c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
1254f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
1255f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
1256f39d1f56SLois Curfman McInnes 
1257ff756334SLois Curfman McInnes    Notes:
1258ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
1259ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
1260ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
1261e51e0e81SBarry Smith 
126295452b02SPatrick Sanan    Fortran Notes:
126395452b02SPatrick Sanan     To use this from Fortran with a ctx you must write an interface definition for this
1264daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
1265daf670e6SBarry Smith     in as the ctx argument.
1266f38a8d56SBarry Smith 
1267f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
1268f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
1269645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
1270645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
1271645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
1272645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
1273645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
1274645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
1275645985a0SLois Curfman McInnes    For example,
1276f39d1f56SLois Curfman McInnes 
1277f39d1f56SLois Curfman McInnes $
1278f39d1f56SLois Curfman McInnes $     Vec x, y
12795bd1e576SStefano Zampini $     extern PetscErrorCode mult(Mat,Vec,Vec);
1280645985a0SLois Curfman McInnes $     Mat A
1281f39d1f56SLois Curfman McInnes $
1282c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
1283c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
1284f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
1285c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
1286c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
1287c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
1288645985a0SLois Curfman McInnes $     MatMult(A,x,y);
128945960306SStefano Zampini $     MatDestroy(&A);
129045960306SStefano Zampini $     VecDestroy(&y);
129145960306SStefano Zampini $     VecDestroy(&x);
1292645985a0SLois Curfman McInnes $
1293e51e0e81SBarry Smith 
12949b53d762SBarry Smith 
129545960306SStefano Zampini    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), MatScale(), MatZeroRows() and MatZeroRowsColumns() internally, so these
12969b53d762SBarry Smith    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.
12979b53d762SBarry Smith 
12989b53d762SBarry Smith 
12990c0fd78eSBarry Smith     For rectangular matrices do all the scalings and shifts make sense?
13000c0fd78eSBarry Smith 
130195452b02SPatrick Sanan     Developers Notes:
130295452b02SPatrick Sanan     Regarding shifting and scaling. The general form is
13030c0fd78eSBarry Smith 
13040c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
13050c0fd78eSBarry Smith 
13060c0fd78eSBarry Smith       The order you apply the operations is important. For example if you have a dshift then
13070c0fd78eSBarry Smith       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
13080c0fd78eSBarry Smith       you get s*vscale*A + diag(shift)
13090c0fd78eSBarry Smith 
13100c0fd78eSBarry Smith           A is the user provided function.
13110c0fd78eSBarry Smith 
1312ad2f5c3fSBarry Smith    KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for
1313ad2f5c3fSBarry Smith    for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger
1314ad2f5c3fSBarry Smith    an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat);
1315ad2f5c3fSBarry Smith    each time the MATSHELL matrix has changed.
1316ad2f5c3fSBarry Smith 
1317ad2f5c3fSBarry Smith    Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided
1318ad2f5c3fSBarry Smith    with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale().
1319ad2f5c3fSBarry Smith 
13200c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
1321e51e0e81SBarry Smith @*/
13227087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
1323e51e0e81SBarry Smith {
1324dfbe8321SBarry Smith   PetscErrorCode ierr;
1325ed3cc1f0SBarry Smith 
13263a40ed3dSBarry Smith   PetscFunctionBegin;
1327f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1328f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1329273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
1330273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
13310fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
1332273d9f13SBarry Smith   PetscFunctionReturn(0);
1333c7fcc2eaSBarry Smith }
1334c7fcc2eaSBarry Smith 
1335789d8953SBarry Smith 
1336c6866cfdSSatish Balay /*@
1337273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
1338c7fcc2eaSBarry Smith 
13393f9fe445SBarry Smith    Logically Collective on Mat
1340c7fcc2eaSBarry Smith 
1341273d9f13SBarry Smith     Input Parameters:
1342273d9f13SBarry Smith +   mat - the shell matrix
1343273d9f13SBarry Smith -   ctx - the context
1344273d9f13SBarry Smith 
1345273d9f13SBarry Smith    Level: advanced
1346273d9f13SBarry Smith 
134795452b02SPatrick Sanan    Fortran Notes:
134895452b02SPatrick Sanan     To use this from Fortran you must write a Fortran interface definition for this
1349daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1350273d9f13SBarry Smith 
1351273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
13520bc0a719SBarry Smith @*/
13537087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
1354273d9f13SBarry Smith {
1355dfbe8321SBarry Smith   PetscErrorCode ierr;
1356273d9f13SBarry Smith 
1357273d9f13SBarry Smith   PetscFunctionBegin;
13580700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1359*ef5c7bd2SStefano Zampini   ierr = PetscTryMethod(mat,"MatShellSetContext_C",(Mat,void*),(mat,ctx));CHKERRQ(ierr);
13603a40ed3dSBarry Smith   PetscFunctionReturn(0);
1361e51e0e81SBarry Smith }
1362e51e0e81SBarry Smith 
1363db77db73SJeremy L Thompson /*@C
1364db77db73SJeremy L Thompson  MatShellSetVecType - Sets the type of Vec returned by MatCreateVecs()
1365db77db73SJeremy L Thompson 
1366db77db73SJeremy L Thompson  Logically collective
1367db77db73SJeremy L Thompson 
1368db77db73SJeremy L Thompson     Input Parameters:
1369db77db73SJeremy L Thompson +   mat   - the shell matrix
1370db77db73SJeremy L Thompson -   vtype - type to use for creating vectors
1371db77db73SJeremy L Thompson 
1372db77db73SJeremy L Thompson  Notes:
1373db77db73SJeremy L Thompson 
1374db77db73SJeremy L Thompson  Level: advanced
1375db77db73SJeremy L Thompson 
1376db77db73SJeremy L Thompson .seealso: MatCreateVecs()
1377db77db73SJeremy L Thompson @*/
1378db77db73SJeremy L Thompson PetscErrorCode  MatShellSetVecType(Mat mat,VecType vtype)
1379db77db73SJeremy L Thompson {
1380db77db73SJeremy L Thompson   PetscErrorCode ierr;
1381db77db73SJeremy L Thompson 
1382db77db73SJeremy L Thompson   PetscFunctionBegin;
1383db77db73SJeremy L Thompson   ierr = PetscTryMethod(mat,"MatShellSetVecType_C",(Mat,VecType),(mat,vtype));CHKERRQ(ierr);
1384db77db73SJeremy L Thompson   PetscFunctionReturn(0);
1385db77db73SJeremy L Thompson }
1386db77db73SJeremy L Thompson 
13870c0fd78eSBarry Smith /*@
13880c0fd78eSBarry Smith     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
13890c0fd78eSBarry Smith           after MatCreateShell()
13900c0fd78eSBarry Smith 
13910c0fd78eSBarry Smith    Logically Collective on Mat
13920c0fd78eSBarry Smith 
13930c0fd78eSBarry Smith     Input Parameter:
13940c0fd78eSBarry Smith .   mat - the shell matrix
13950c0fd78eSBarry Smith 
13960c0fd78eSBarry Smith   Level: advanced
13970c0fd78eSBarry Smith 
13980c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
13990c0fd78eSBarry Smith @*/
14000c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A)
14010c0fd78eSBarry Smith {
14020c0fd78eSBarry Smith   PetscErrorCode ierr;
14030c0fd78eSBarry Smith 
14040c0fd78eSBarry Smith   PetscFunctionBegin;
14050c0fd78eSBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1406*ef5c7bd2SStefano Zampini   ierr = PetscTryMethod(A,"MatShellSetManageScalingShifts_C",(Mat),(A));CHKERRQ(ierr);
14070c0fd78eSBarry Smith   PetscFunctionReturn(0);
14080c0fd78eSBarry Smith }
14090c0fd78eSBarry Smith 
1410c16cb8f2SBarry Smith /*@C
1411f3b1f45cSBarry Smith     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
1412f3b1f45cSBarry Smith 
1413f3b1f45cSBarry Smith    Logically Collective on Mat
1414f3b1f45cSBarry Smith 
1415f3b1f45cSBarry Smith     Input Parameters:
1416f3b1f45cSBarry Smith +   mat - the shell matrix
1417f3b1f45cSBarry Smith .   f - the function
1418f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
1419f3b1f45cSBarry Smith -   ctx - an optional context for the function
1420f3b1f45cSBarry Smith 
1421f3b1f45cSBarry Smith    Output Parameter:
1422f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
1423f3b1f45cSBarry Smith 
1424f3b1f45cSBarry Smith    Options Database:
1425f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
1426f3b1f45cSBarry Smith 
1427f3b1f45cSBarry Smith    Level: advanced
1428f3b1f45cSBarry Smith 
142995452b02SPatrick Sanan    Fortran Notes:
143095452b02SPatrick Sanan     Not supported from Fortran
1431f3b1f45cSBarry Smith 
1432f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
1433f3b1f45cSBarry Smith @*/
1434f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
1435f3b1f45cSBarry Smith {
1436f3b1f45cSBarry Smith   PetscErrorCode ierr;
1437f3b1f45cSBarry Smith   PetscInt       m,n;
1438f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
1439f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
144074e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
1441f3b1f45cSBarry Smith 
1442f3b1f45cSBarry Smith   PetscFunctionBegin;
1443f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1444f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr);
1445f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
1446f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr);
1447f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
1448f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
1449f3b1f45cSBarry Smith 
14500bacdadaSStefano Zampini   ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr);
14510bacdadaSStefano Zampini   ierr = MatComputeOperator(mat,MATAIJ,&Dmat);CHKERRQ(ierr);
1452f3b1f45cSBarry Smith 
1453f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
1454f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
1455f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
1456f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
1457f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
1458f3b1f45cSBarry Smith     flag = PETSC_FALSE;
1459f3b1f45cSBarry Smith     if (v) {
1460fc7aafd1SBarry Smith       ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm));CHKERRQ(ierr);
1461f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
1462f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
1463f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
1464f3b1f45cSBarry Smith     }
1465f3b1f45cSBarry Smith   } else if (v) {
1466fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
1467f3b1f45cSBarry Smith   }
1468f3b1f45cSBarry Smith   if (flg) *flg = flag;
1469f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
1470f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
1471f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
1472f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
1473f3b1f45cSBarry Smith   PetscFunctionReturn(0);
1474f3b1f45cSBarry Smith }
1475f3b1f45cSBarry Smith 
1476f3b1f45cSBarry Smith /*@C
1477f3b1f45cSBarry Smith     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.
1478f3b1f45cSBarry Smith 
1479f3b1f45cSBarry Smith    Logically Collective on Mat
1480f3b1f45cSBarry Smith 
1481f3b1f45cSBarry Smith     Input Parameters:
1482f3b1f45cSBarry Smith +   mat - the shell matrix
1483f3b1f45cSBarry Smith .   f - the function
1484f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
1485f3b1f45cSBarry Smith -   ctx - an optional context for the function
1486f3b1f45cSBarry Smith 
1487f3b1f45cSBarry Smith    Output Parameter:
1488f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
1489f3b1f45cSBarry Smith 
1490f3b1f45cSBarry Smith    Options Database:
1491f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
1492f3b1f45cSBarry Smith 
1493f3b1f45cSBarry Smith    Level: advanced
1494f3b1f45cSBarry Smith 
149595452b02SPatrick Sanan    Fortran Notes:
149695452b02SPatrick Sanan     Not supported from Fortran
1497f3b1f45cSBarry Smith 
1498f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
1499f3b1f45cSBarry Smith @*/
1500f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
1501f3b1f45cSBarry Smith {
1502f3b1f45cSBarry Smith   PetscErrorCode ierr;
1503f3b1f45cSBarry Smith   Vec            x,y,z;
1504f3b1f45cSBarry Smith   PetscInt       m,n,M,N;
1505f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
1506f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
150774e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
1508f3b1f45cSBarry Smith 
1509f3b1f45cSBarry Smith   PetscFunctionBegin;
1510f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1511f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
1512f3b1f45cSBarry Smith   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
1513f3b1f45cSBarry Smith   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
1514f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
1515f3b1f45cSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
1516f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
1517f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
1518f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
15190bacdadaSStefano Zampini   ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr);
1520f3b1f45cSBarry Smith   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
15210bacdadaSStefano Zampini   ierr = MatComputeOperatorTranspose(mat,MATAIJ,&Dmat);CHKERRQ(ierr);
1522f3b1f45cSBarry Smith 
1523f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
1524f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
1525f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
1526f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
1527f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
1528f3b1f45cSBarry Smith     flag = PETSC_FALSE;
1529f3b1f45cSBarry Smith     if (v) {
1530fc7aafd1SBarry Smith       ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm));CHKERRQ(ierr);
1531f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
1532f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
1533f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
1534f3b1f45cSBarry Smith     }
1535f3b1f45cSBarry Smith   } else if (v) {
1536fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
1537f3b1f45cSBarry Smith   }
1538f3b1f45cSBarry Smith   if (flg) *flg = flag;
1539f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
1540f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
1541f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
1542f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
1543f3b1f45cSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
1544f3b1f45cSBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
1545f3b1f45cSBarry Smith   ierr = VecDestroy(&z);CHKERRQ(ierr);
1546f3b1f45cSBarry Smith   PetscFunctionReturn(0);
1547f3b1f45cSBarry Smith }
1548f3b1f45cSBarry Smith 
1549f3b1f45cSBarry Smith /*@C
15500c0fd78eSBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
1551e51e0e81SBarry Smith 
15523f9fe445SBarry Smith    Logically Collective on Mat
1553fee21e36SBarry Smith 
1554c7fcc2eaSBarry Smith     Input Parameters:
1555c7fcc2eaSBarry Smith +   mat - the shell matrix
1556c7fcc2eaSBarry Smith .   op - the name of the operation
1557789d8953SBarry Smith -   g - the function that provides the operation.
1558c7fcc2eaSBarry Smith 
155915091d37SBarry Smith    Level: advanced
156015091d37SBarry Smith 
1561fae171e0SBarry Smith     Usage:
1562e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
1563f39d1f56SLois Curfman McInnes $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
1564c134de8dSSatish Balay $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
15650b627109SLois Curfman McInnes 
1566a62d957aSLois Curfman McInnes     Notes:
1567e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
15681c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
1569a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
15701c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
1571a62d957aSLois Curfman McInnes 
1572a4d85fd7SBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
1573deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
1574deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
1575deebb3c3SLois Curfman McInnes     routines, e.g.,
1576a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1577a62d957aSLois Curfman McInnes 
15784aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
15794aa34b0aSBarry Smith     nonzero on failure.
15804aa34b0aSBarry Smith 
1581a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
1582a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
1583a62d957aSLois Curfman McInnes     set by MatCreateShell().
1584a62d957aSLois Curfman McInnes 
158595452b02SPatrick Sanan     Fortran Notes:
158695452b02SPatrick Sanan     For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
1587501d9185SBarry Smith        generate a matrix. See src/mat/examples/tests/ex120f.F
1588501d9185SBarry Smith 
15890c0fd78eSBarry Smith     Use MatSetOperation() to set an operation for any matrix type
15900c0fd78eSBarry Smith 
15910c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts()
1592e51e0e81SBarry Smith @*/
1593789d8953SBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*g)(void))
1594e51e0e81SBarry Smith {
1595976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
1596273d9f13SBarry Smith 
15973a40ed3dSBarry Smith   PetscFunctionBegin;
15980700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1599*ef5c7bd2SStefano Zampini   ierr = PetscTryMethod(mat,"MatShellSetOperation_C",(Mat,MatOperation,void (*)(void)),(mat,op,g));CHKERRQ(ierr);
16003a40ed3dSBarry Smith   PetscFunctionReturn(0);
1601e51e0e81SBarry Smith }
1602f0479e8cSBarry Smith 
1603d4bb536fSBarry Smith /*@C
1604d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
1605d4bb536fSBarry Smith 
1606c7fcc2eaSBarry Smith     Not Collective
1607c7fcc2eaSBarry Smith 
1608d4bb536fSBarry Smith     Input Parameters:
1609c7fcc2eaSBarry Smith +   mat - the shell matrix
1610c7fcc2eaSBarry Smith -   op - the name of the operation
1611d4bb536fSBarry Smith 
1612d4bb536fSBarry Smith     Output Parameter:
1613789d8953SBarry Smith .   g - the function that provides the operation.
1614d4bb536fSBarry Smith 
161515091d37SBarry Smith     Level: advanced
161615091d37SBarry Smith 
1617d4bb536fSBarry Smith     Notes:
1618e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
1619d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
1620d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
1621d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
1622d4bb536fSBarry Smith 
1623d4bb536fSBarry Smith     All user-provided functions have the same calling
1624d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
1625d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
1626d4bb536fSBarry Smith     routines, e.g.,
1627d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1628d4bb536fSBarry Smith 
1629d4bb536fSBarry Smith     Within each user-defined routine, the user should call
1630d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
1631d4bb536fSBarry Smith     set by MatCreateShell().
1632d4bb536fSBarry Smith 
1633ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1634d4bb536fSBarry Smith @*/
1635789d8953SBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**g)(void))
1636d4bb536fSBarry Smith {
1637976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
1638273d9f13SBarry Smith 
16393a40ed3dSBarry Smith   PetscFunctionBegin;
16400700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1641789d8953SBarry Smith   ierr = PetscUseMethod(mat,"MatShellGetOperation_C",(Mat,MatOperation,void (**)(void)),(mat,op,g));CHKERRQ(ierr);
16423a40ed3dSBarry Smith   PetscFunctionReturn(0);
1643d4bb536fSBarry Smith }
1644