xref: /petsc/src/mat/impls/shell/shell.c (revision 843d480f9cddad1d931d742d77ea89a67ca9f317)
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 
19b77ba244SStefano Zampini struct _n_MatShellMatFunctionList {
20b77ba244SStefano Zampini   PetscErrorCode  (*symbolic)(Mat,Mat,Mat,void**);
21b77ba244SStefano Zampini   PetscErrorCode  (*numeric)(Mat,Mat,Mat,void*);
22b77ba244SStefano Zampini   PetscErrorCode  (*destroy)(void*);
23b77ba244SStefano Zampini   MatProductType  ptype;
24b77ba244SStefano Zampini   char            *composedname;  /* string to identify routine with double dispatch */
25b77ba244SStefano Zampini   char            *resultname; /* result matrix type */
26b77ba244SStefano Zampini 
27b77ba244SStefano Zampini   struct _n_MatShellMatFunctionList *next;
28b77ba244SStefano Zampini };
29b77ba244SStefano Zampini typedef struct _n_MatShellMatFunctionList *MatShellMatFunctionList;
30b77ba244SStefano Zampini 
3128f357e3SAlex Fikl typedef struct {
3228f357e3SAlex Fikl   struct _MatShellOps ops[1];
332205254eSKarl Rupp 
34b77ba244SStefano Zampini   /* The user will manage the scaling and shifts for the MATSHELL, not the default */
35b77ba244SStefano Zampini   PetscBool managescalingshifts;
36b77ba244SStefano Zampini 
37b77ba244SStefano Zampini   /* support for MatScale, MatShift and MatMultAdd */
38ef66eb69SBarry Smith   PetscScalar vscale,vshift;
398fe8eb27SJed Brown   Vec         dshift;
408fe8eb27SJed Brown   Vec         left,right;
418fe8eb27SJed Brown   Vec         left_work,right_work;
425edf6516SJed Brown   Vec         left_add_work,right_add_work;
43b77ba244SStefano Zampini 
44ef5c7bd2SStefano Zampini   /* support for MatAXPY */
459f137db4SBarry Smith   Mat              axpy;
469f137db4SBarry Smith   PetscScalar      axpy_vscale;
47b77ba244SStefano Zampini   Vec              axpy_left,axpy_right;
48ef5c7bd2SStefano Zampini   PetscObjectState axpy_state;
49b77ba244SStefano Zampini 
5045960306SStefano Zampini   /* support for ZeroRows/Columns operations */
5145960306SStefano Zampini   IS         zrows;
5245960306SStefano Zampini   IS         zcols;
5345960306SStefano Zampini   Vec        zvals;
5445960306SStefano Zampini   Vec        zvals_w;
5545960306SStefano Zampini   VecScatter zvals_sct_r;
5645960306SStefano Zampini   VecScatter zvals_sct_c;
57b77ba244SStefano Zampini 
58b77ba244SStefano Zampini   /* MatMat operations */
59b77ba244SStefano Zampini   MatShellMatFunctionList matmat;
60b77ba244SStefano Zampini 
61b77ba244SStefano Zampini   /* user defined context */
6220563c6bSBarry Smith   void *ctx;
6388cf3e7dSBarry Smith } Mat_Shell;
640c0fd78eSBarry Smith 
6545960306SStefano Zampini 
6645960306SStefano Zampini /*
6745960306SStefano Zampini      Store and scale values on zeroed rows
6845960306SStefano Zampini      xx = [x_1, 0], 0 on zeroed columns
6945960306SStefano Zampini */
7045960306SStefano Zampini static PetscErrorCode MatShellPreZeroRight(Mat A,Vec x,Vec *xx)
7145960306SStefano Zampini {
7245960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
7345960306SStefano Zampini   PetscErrorCode ierr;
7445960306SStefano Zampini 
7545960306SStefano Zampini   PetscFunctionBegin;
7645960306SStefano Zampini   *xx = x;
7745960306SStefano Zampini   if (shell->zrows) {
7845960306SStefano Zampini     ierr = VecSet(shell->zvals_w,0.0);CHKERRQ(ierr);
7945960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8045960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8145960306SStefano Zampini     ierr = VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);CHKERRQ(ierr);
8245960306SStefano Zampini   }
8345960306SStefano Zampini   if (shell->zcols) {
8445960306SStefano Zampini     if (!shell->right_work) {
8545960306SStefano Zampini       ierr = MatCreateVecs(A,&shell->right_work,NULL);CHKERRQ(ierr);
8645960306SStefano Zampini     }
8745960306SStefano Zampini     ierr = VecCopy(x,shell->right_work);CHKERRQ(ierr);
8845960306SStefano Zampini     ierr = VecISSet(shell->right_work,shell->zcols,0.0);CHKERRQ(ierr);
8945960306SStefano Zampini     *xx  = shell->right_work;
9045960306SStefano Zampini   }
9145960306SStefano Zampini   PetscFunctionReturn(0);
9245960306SStefano Zampini }
9345960306SStefano Zampini 
9445960306SStefano Zampini /* Insert properly diagonally scaled values stored in MatShellPreZeroRight */
9545960306SStefano Zampini static PetscErrorCode MatShellPostZeroLeft(Mat A,Vec x)
9645960306SStefano Zampini {
9745960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
9845960306SStefano Zampini   PetscErrorCode ierr;
9945960306SStefano Zampini 
10045960306SStefano Zampini   PetscFunctionBegin;
10145960306SStefano Zampini   if (shell->zrows) {
10245960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10345960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10445960306SStefano Zampini   }
10545960306SStefano Zampini   PetscFunctionReturn(0);
10645960306SStefano Zampini }
10745960306SStefano Zampini 
10845960306SStefano Zampini /*
10945960306SStefano Zampini      Store and scale values on zeroed rows
11045960306SStefano Zampini      xx = [x_1, 0], 0 on zeroed rows
11145960306SStefano Zampini */
11245960306SStefano Zampini static PetscErrorCode MatShellPreZeroLeft(Mat A,Vec x,Vec *xx)
11345960306SStefano Zampini {
11445960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
11545960306SStefano Zampini   PetscErrorCode ierr;
11645960306SStefano Zampini 
11745960306SStefano Zampini   PetscFunctionBegin;
11845960306SStefano Zampini   *xx = NULL;
11945960306SStefano Zampini   if (!shell->zrows) {
12045960306SStefano Zampini     *xx = x;
12145960306SStefano Zampini   } else {
12245960306SStefano Zampini     if (!shell->left_work) {
12345960306SStefano Zampini       ierr = MatCreateVecs(A,NULL,&shell->left_work);CHKERRQ(ierr);
12445960306SStefano Zampini     }
12545960306SStefano Zampini     ierr = VecCopy(x,shell->left_work);CHKERRQ(ierr);
12645960306SStefano Zampini     ierr = VecSet(shell->zvals_w,0.0);CHKERRQ(ierr);
12745960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12845960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12945960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13045960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13145960306SStefano Zampini     ierr = VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);CHKERRQ(ierr);
13245960306SStefano Zampini     *xx  = shell->left_work;
13345960306SStefano Zampini   }
13445960306SStefano Zampini   PetscFunctionReturn(0);
13545960306SStefano Zampini }
13645960306SStefano Zampini 
13745960306SStefano Zampini /* Zero zero-columns contributions, sum contributions from properly scaled values stored in MatShellPreZeroLeft */
13845960306SStefano Zampini static PetscErrorCode MatShellPostZeroRight(Mat A,Vec x)
13945960306SStefano Zampini {
14045960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
14145960306SStefano Zampini   PetscErrorCode ierr;
14245960306SStefano Zampini 
14345960306SStefano Zampini   PetscFunctionBegin;
14445960306SStefano Zampini   if (shell->zcols) {
14545960306SStefano Zampini     ierr = VecISSet(x,shell->zcols,0.0);CHKERRQ(ierr);
14645960306SStefano Zampini   }
14745960306SStefano Zampini   if (shell->zrows) {
14845960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14945960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15045960306SStefano Zampini   }
15145960306SStefano Zampini   PetscFunctionReturn(0);
15245960306SStefano Zampini }
15345960306SStefano Zampini 
1548fe8eb27SJed Brown /*
1550c0fd78eSBarry Smith       xx = diag(left)*x
1568fe8eb27SJed Brown */
1578fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
1588fe8eb27SJed Brown {
1598fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1608fe8eb27SJed Brown   PetscErrorCode ierr;
1618fe8eb27SJed Brown 
1628fe8eb27SJed Brown   PetscFunctionBegin;
1630298fd71SBarry Smith   *xx = NULL;
1648fe8eb27SJed Brown   if (!shell->left) {
1658fe8eb27SJed Brown     *xx = x;
1668fe8eb27SJed Brown   } else {
1678fe8eb27SJed Brown     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
1688fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
1698fe8eb27SJed Brown     *xx  = shell->left_work;
1708fe8eb27SJed Brown   }
1718fe8eb27SJed Brown   PetscFunctionReturn(0);
1728fe8eb27SJed Brown }
1738fe8eb27SJed Brown 
1740c0fd78eSBarry Smith /*
1750c0fd78eSBarry Smith      xx = diag(right)*x
1760c0fd78eSBarry Smith */
1778fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
1788fe8eb27SJed Brown {
1798fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
1808fe8eb27SJed Brown   PetscErrorCode ierr;
1818fe8eb27SJed Brown 
1828fe8eb27SJed Brown   PetscFunctionBegin;
1830298fd71SBarry Smith   *xx = NULL;
1848fe8eb27SJed Brown   if (!shell->right) {
1858fe8eb27SJed Brown     *xx = x;
1868fe8eb27SJed Brown   } else {
1878fe8eb27SJed Brown     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
1888fe8eb27SJed Brown     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
1898fe8eb27SJed Brown     *xx  = shell->right_work;
1908fe8eb27SJed Brown   }
1918fe8eb27SJed Brown   PetscFunctionReturn(0);
1928fe8eb27SJed Brown }
1938fe8eb27SJed Brown 
1940c0fd78eSBarry Smith /*
1950c0fd78eSBarry Smith     x = diag(left)*x
1960c0fd78eSBarry Smith */
1978fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
1988fe8eb27SJed Brown {
1998fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
2008fe8eb27SJed Brown   PetscErrorCode ierr;
2018fe8eb27SJed Brown 
2028fe8eb27SJed Brown   PetscFunctionBegin;
2038fe8eb27SJed Brown   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
2048fe8eb27SJed Brown   PetscFunctionReturn(0);
2058fe8eb27SJed Brown }
2068fe8eb27SJed Brown 
2070c0fd78eSBarry Smith /*
2080c0fd78eSBarry Smith     x = diag(right)*x
2090c0fd78eSBarry Smith */
2108fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
2118fe8eb27SJed Brown {
2128fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
2138fe8eb27SJed Brown   PetscErrorCode ierr;
2148fe8eb27SJed Brown 
2158fe8eb27SJed Brown   PetscFunctionBegin;
2168fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
2178fe8eb27SJed Brown   PetscFunctionReturn(0);
2188fe8eb27SJed Brown }
2198fe8eb27SJed Brown 
2200c0fd78eSBarry Smith /*
2210c0fd78eSBarry Smith          Y = vscale*Y + diag(dshift)*X + vshift*X
2220c0fd78eSBarry Smith 
2230c0fd78eSBarry Smith          On input Y already contains A*x
2240c0fd78eSBarry Smith */
2258fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
2268fe8eb27SJed Brown {
2278fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
2288fe8eb27SJed Brown   PetscErrorCode ierr;
2298fe8eb27SJed Brown 
2308fe8eb27SJed Brown   PetscFunctionBegin;
2318fe8eb27SJed Brown   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
2328fe8eb27SJed Brown     PetscInt          i,m;
2338fe8eb27SJed Brown     const PetscScalar *x,*d;
2348fe8eb27SJed Brown     PetscScalar       *y;
2358fe8eb27SJed Brown     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
2368fe8eb27SJed Brown     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
2378fe8eb27SJed Brown     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
2388fe8eb27SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
2398fe8eb27SJed Brown     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
2408fe8eb27SJed Brown     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
2418fe8eb27SJed Brown     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
2428fe8eb27SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
2430c0fd78eSBarry Smith   } else {
244027c5fb5SJed Brown     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
2458fe8eb27SJed Brown   }
246d4c7de66SBarry Smith   if (shell->vshift != 0.0) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */
2478fe8eb27SJed Brown   PetscFunctionReturn(0);
2488fe8eb27SJed Brown }
249e51e0e81SBarry Smith 
250789d8953SBarry Smith PetscErrorCode MatShellGetContext_Shell(Mat mat,void *ctx)
251789d8953SBarry Smith {
252789d8953SBarry Smith   PetscFunctionBegin;
253789d8953SBarry Smith   *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
254789d8953SBarry Smith   PetscFunctionReturn(0);
255789d8953SBarry Smith }
256789d8953SBarry Smith 
2579d225801SJed Brown /*@
258a62d957aSLois Curfman McInnes     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
259b4fd4287SBarry Smith 
260c7fcc2eaSBarry Smith     Not Collective
261c7fcc2eaSBarry Smith 
262b4fd4287SBarry Smith     Input Parameter:
263b4fd4287SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
264b4fd4287SBarry Smith 
265b4fd4287SBarry Smith     Output Parameter:
266b4fd4287SBarry Smith .   ctx - the user provided context
267b4fd4287SBarry Smith 
26815091d37SBarry Smith     Level: advanced
26915091d37SBarry Smith 
27095452b02SPatrick Sanan    Fortran Notes:
27195452b02SPatrick Sanan     To use this from Fortran you must write a Fortran interface definition for this
272daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
273a62d957aSLois Curfman McInnes 
274ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
2750bc0a719SBarry Smith @*/
2768fe8eb27SJed Brown PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
277b4fd4287SBarry Smith {
278dfbe8321SBarry Smith   PetscErrorCode ierr;
279273d9f13SBarry Smith 
2803a40ed3dSBarry Smith   PetscFunctionBegin;
2810700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2824482741eSBarry Smith   PetscValidPointer(ctx,2);
283789d8953SBarry Smith   ierr = PetscUseMethod(mat,"MatShellGetContext_C",(Mat,void*),(mat,ctx));CHKERRQ(ierr);
2843a40ed3dSBarry Smith   PetscFunctionReturn(0);
285b4fd4287SBarry Smith }
286b4fd4287SBarry Smith 
28745960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Local_Shell(Mat mat,PetscInt nr,PetscInt rows[],PetscInt nc,PetscInt cols[],PetscScalar diag,PetscBool rc)
28845960306SStefano Zampini {
28945960306SStefano Zampini   PetscErrorCode ierr;
29045960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)mat->data;
29145960306SStefano Zampini   Vec            x = NULL,b = NULL;
29245960306SStefano Zampini   IS             is1, is2;
29345960306SStefano Zampini   const PetscInt *ridxs;
29445960306SStefano Zampini   PetscInt       *idxs,*gidxs;
29545960306SStefano Zampini   PetscInt       cum,rst,cst,i;
29645960306SStefano Zampini 
29745960306SStefano Zampini   PetscFunctionBegin;
29845960306SStefano Zampini   if (!shell->zvals) {
29945960306SStefano Zampini     ierr = MatCreateVecs(mat,NULL,&shell->zvals);CHKERRQ(ierr);
30045960306SStefano Zampini   }
30145960306SStefano Zampini   if (!shell->zvals_w) {
30245960306SStefano Zampini     ierr = VecDuplicate(shell->zvals,&shell->zvals_w);CHKERRQ(ierr);
30345960306SStefano Zampini   }
30445960306SStefano Zampini   ierr = MatGetOwnershipRange(mat,&rst,NULL);CHKERRQ(ierr);
30545960306SStefano Zampini   ierr = MatGetOwnershipRangeColumn(mat,&cst,NULL);CHKERRQ(ierr);
30645960306SStefano Zampini 
30745960306SStefano Zampini   /* Expand/create index set of zeroed rows */
30845960306SStefano Zampini   ierr = PetscMalloc1(nr,&idxs);CHKERRQ(ierr);
30945960306SStefano Zampini   for (i = 0; i < nr; i++) idxs[i] = rows[i] + rst;
31045960306SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,nr,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr);
31145960306SStefano Zampini   ierr = ISSort(is1);CHKERRQ(ierr);
31245960306SStefano Zampini   ierr = VecISSet(shell->zvals,is1,diag);CHKERRQ(ierr);
31345960306SStefano Zampini   if (shell->zrows) {
31445960306SStefano Zampini     ierr = ISSum(shell->zrows,is1,&is2);CHKERRQ(ierr);
31545960306SStefano Zampini     ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr);
31645960306SStefano Zampini     ierr = ISDestroy(&is1);CHKERRQ(ierr);
31745960306SStefano Zampini     shell->zrows = is2;
31845960306SStefano Zampini   } else shell->zrows = is1;
31945960306SStefano Zampini 
32045960306SStefano Zampini   /* Create scatters for diagonal values communications */
32145960306SStefano Zampini   ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr);
32245960306SStefano Zampini   ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr);
32345960306SStefano Zampini 
32445960306SStefano Zampini   /* row scatter: from/to left vector */
32545960306SStefano Zampini   ierr = MatCreateVecs(mat,&x,&b);CHKERRQ(ierr);
32645960306SStefano Zampini   ierr = VecScatterCreate(b,shell->zrows,shell->zvals_w,shell->zrows,&shell->zvals_sct_r);CHKERRQ(ierr);
32745960306SStefano Zampini 
32845960306SStefano Zampini   /* col scatter: from right vector to left vector */
32945960306SStefano Zampini   ierr = ISGetIndices(shell->zrows,&ridxs);CHKERRQ(ierr);
33045960306SStefano Zampini   ierr = ISGetLocalSize(shell->zrows,&nr);CHKERRQ(ierr);
33145960306SStefano Zampini   ierr = PetscMalloc1(nr,&gidxs);CHKERRQ(ierr);
33245960306SStefano Zampini   for (i = 0, cum  = 0; i < nr; i++) {
33345960306SStefano Zampini     if (ridxs[i] >= mat->cmap->N) continue;
33445960306SStefano Zampini     gidxs[cum] = ridxs[i];
33545960306SStefano Zampini     cum++;
33645960306SStefano Zampini   }
33745960306SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)mat),cum,gidxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr);
33845960306SStefano Zampini   ierr = VecScatterCreate(x,is1,shell->zvals_w,is1,&shell->zvals_sct_c);CHKERRQ(ierr);
33945960306SStefano Zampini   ierr = ISDestroy(&is1);CHKERRQ(ierr);
34045960306SStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
34145960306SStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
34245960306SStefano Zampini 
34345960306SStefano Zampini   /* Expand/create index set of zeroed columns */
34445960306SStefano Zampini   if (rc) {
34545960306SStefano Zampini     ierr = PetscMalloc1(nc,&idxs);CHKERRQ(ierr);
34645960306SStefano Zampini     for (i = 0; i < nc; i++) idxs[i] = cols[i] + cst;
34745960306SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nc,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr);
34845960306SStefano Zampini     ierr = ISSort(is1);CHKERRQ(ierr);
34945960306SStefano Zampini     if (shell->zcols) {
35045960306SStefano Zampini       ierr = ISSum(shell->zcols,is1,&is2);CHKERRQ(ierr);
35145960306SStefano Zampini       ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr);
35245960306SStefano Zampini       ierr = ISDestroy(&is1);CHKERRQ(ierr);
35345960306SStefano Zampini       shell->zcols = is2;
35445960306SStefano Zampini     } else shell->zcols = is1;
35545960306SStefano Zampini   }
35645960306SStefano Zampini   PetscFunctionReturn(0);
35745960306SStefano Zampini }
35845960306SStefano Zampini 
35945960306SStefano Zampini static PetscErrorCode MatZeroRows_Shell(Mat mat,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
36045960306SStefano Zampini {
361ef5c7bd2SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)mat->data;
36245960306SStefano Zampini   PetscInt       nr, *lrows;
36345960306SStefano Zampini   PetscErrorCode ierr;
36445960306SStefano Zampini 
36545960306SStefano Zampini   PetscFunctionBegin;
36645960306SStefano Zampini   if (x && b) {
36745960306SStefano Zampini     Vec          xt;
36845960306SStefano Zampini     PetscScalar *vals;
36945960306SStefano Zampini     PetscInt    *gcols,i,st,nl,nc;
37045960306SStefano Zampini 
37145960306SStefano Zampini     ierr = PetscMalloc1(n,&gcols);CHKERRQ(ierr);
37245960306SStefano Zampini     for (i = 0, nc = 0; i < n; i++) if (rows[i] < mat->cmap->N) gcols[nc++] = rows[i];
37345960306SStefano Zampini 
37445960306SStefano Zampini     ierr = MatCreateVecs(mat,&xt,NULL);CHKERRQ(ierr);
37545960306SStefano Zampini     ierr = VecCopy(x,xt);CHKERRQ(ierr);
37645960306SStefano Zampini     ierr = PetscCalloc1(nc,&vals);CHKERRQ(ierr);
37745960306SStefano Zampini     ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */
37845960306SStefano Zampini     ierr = PetscFree(vals);CHKERRQ(ierr);
37945960306SStefano Zampini     ierr = VecAssemblyBegin(xt);CHKERRQ(ierr);
38045960306SStefano Zampini     ierr = VecAssemblyEnd(xt);CHKERRQ(ierr);
38145960306SStefano Zampini     ierr = VecAYPX(xt,-1.0,x);CHKERRQ(ierr);                           /* xt = [0, x2] */
38245960306SStefano Zampini 
38345960306SStefano Zampini     ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr);
38445960306SStefano Zampini     ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr);
38545960306SStefano Zampini     ierr = VecGetArray(xt,&vals);CHKERRQ(ierr);
38645960306SStefano Zampini     for (i = 0; i < nl; i++) {
38745960306SStefano Zampini       PetscInt g = i + st;
38845960306SStefano Zampini       if (g > mat->rmap->N) continue;
38945960306SStefano Zampini       if (PetscAbsScalar(vals[i]) == 0.0) continue;
39045960306SStefano Zampini       ierr = VecSetValue(b,g,diag*vals[i],INSERT_VALUES);CHKERRQ(ierr);
39145960306SStefano Zampini     }
39245960306SStefano Zampini     ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr);
39345960306SStefano Zampini     ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
39445960306SStefano Zampini     ierr = VecAssemblyEnd(b);CHKERRQ(ierr);                            /* b  = [b1, x2 * diag] */
39545960306SStefano Zampini     ierr = VecDestroy(&xt);CHKERRQ(ierr);
39645960306SStefano Zampini     ierr = PetscFree(gcols);CHKERRQ(ierr);
39745960306SStefano Zampini   }
39845960306SStefano Zampini   ierr = PetscLayoutMapLocal(mat->rmap,n,rows,&nr,&lrows,NULL);CHKERRQ(ierr);
39945960306SStefano Zampini   ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,0,NULL,diag,PETSC_FALSE);CHKERRQ(ierr);
400ef5c7bd2SStefano Zampini   if (shell->axpy) {
401ef5c7bd2SStefano Zampini     ierr = MatZeroRows(shell->axpy,n,rows,0.0,NULL,NULL);CHKERRQ(ierr);
402ef5c7bd2SStefano Zampini   }
40345960306SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
40445960306SStefano Zampini   PetscFunctionReturn(0);
40545960306SStefano Zampini }
40645960306SStefano Zampini 
40745960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Shell(Mat mat,PetscInt n,const PetscInt rowscols[],PetscScalar diag,Vec x,Vec b)
40845960306SStefano Zampini {
409ef5c7bd2SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)mat->data;
41045960306SStefano Zampini   PetscInt       *lrows, *lcols;
41145960306SStefano Zampini   PetscInt       nr, nc;
41245960306SStefano Zampini   PetscBool      congruent;
41345960306SStefano Zampini   PetscErrorCode ierr;
41445960306SStefano Zampini 
41545960306SStefano Zampini   PetscFunctionBegin;
41645960306SStefano Zampini   if (x && b) {
41745960306SStefano Zampini     Vec          xt, bt;
41845960306SStefano Zampini     PetscScalar *vals;
41945960306SStefano Zampini     PetscInt    *grows,*gcols,i,st,nl;
42045960306SStefano Zampini 
42145960306SStefano Zampini     ierr = PetscMalloc2(n,&grows,n,&gcols);CHKERRQ(ierr);
42245960306SStefano Zampini     for (i = 0, nr = 0; i < n; i++) if (rowscols[i] < mat->rmap->N) grows[nr++] = rowscols[i];
42345960306SStefano Zampini     for (i = 0, nc = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) gcols[nc++] = rowscols[i];
42445960306SStefano Zampini     ierr = PetscCalloc1(n,&vals);CHKERRQ(ierr);
42545960306SStefano Zampini 
42645960306SStefano Zampini     ierr = MatCreateVecs(mat,&xt,&bt);CHKERRQ(ierr);
42745960306SStefano Zampini     ierr = VecCopy(x,xt);CHKERRQ(ierr);
42845960306SStefano Zampini     ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */
42945960306SStefano Zampini     ierr = VecAssemblyBegin(xt);CHKERRQ(ierr);
43045960306SStefano Zampini     ierr = VecAssemblyEnd(xt);CHKERRQ(ierr);
43145960306SStefano Zampini     ierr = VecAXPY(xt,-1.0,x);CHKERRQ(ierr);                           /* xt = [0, -x2] */
43245960306SStefano Zampini     ierr = MatMult(mat,xt,bt);CHKERRQ(ierr);                           /* bt = [-A12*x2,-A22*x2] */
43345960306SStefano Zampini     ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* bt = [-A12*x2,0] */
43445960306SStefano Zampini     ierr = VecAssemblyBegin(bt);CHKERRQ(ierr);
43545960306SStefano Zampini     ierr = VecAssemblyEnd(bt);CHKERRQ(ierr);
43645960306SStefano Zampini     ierr = VecAXPY(b,1.0,bt);CHKERRQ(ierr);                            /* b  = [b1 - A12*x2, b2] */
43745960306SStefano Zampini     ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* b  = [b1 - A12*x2, 0] */
43845960306SStefano Zampini     ierr = VecAssemblyBegin(bt);CHKERRQ(ierr);
43945960306SStefano Zampini     ierr = VecAssemblyEnd(bt);CHKERRQ(ierr);
44045960306SStefano Zampini     ierr = PetscFree(vals);CHKERRQ(ierr);
44145960306SStefano Zampini 
44245960306SStefano Zampini     ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr);
44345960306SStefano Zampini     ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr);
44445960306SStefano Zampini     ierr = VecGetArray(xt,&vals);CHKERRQ(ierr);
44545960306SStefano Zampini     for (i = 0; i < nl; i++) {
44645960306SStefano Zampini       PetscInt g = i + st;
44745960306SStefano Zampini       if (g > mat->rmap->N) continue;
44845960306SStefano Zampini       if (PetscAbsScalar(vals[i]) == 0.0) continue;
44945960306SStefano Zampini       ierr = VecSetValue(b,g,-diag*vals[i],INSERT_VALUES);CHKERRQ(ierr);
45045960306SStefano Zampini     }
45145960306SStefano Zampini     ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr);
45245960306SStefano Zampini     ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
45345960306SStefano Zampini     ierr = VecAssemblyEnd(b);CHKERRQ(ierr);                            /* b  = [b1 - A12*x2, x2 * diag] */
45445960306SStefano Zampini     ierr = VecDestroy(&xt);CHKERRQ(ierr);
45545960306SStefano Zampini     ierr = VecDestroy(&bt);CHKERRQ(ierr);
45645960306SStefano Zampini     ierr = PetscFree2(grows,gcols);CHKERRQ(ierr);
45745960306SStefano Zampini   }
45845960306SStefano Zampini   ierr = PetscLayoutMapLocal(mat->rmap,n,rowscols,&nr,&lrows,NULL);CHKERRQ(ierr);
45945960306SStefano Zampini   ierr = MatHasCongruentLayouts(mat,&congruent);CHKERRQ(ierr);
46045960306SStefano Zampini   if (congruent) {
46145960306SStefano Zampini     nc    = nr;
46245960306SStefano Zampini     lcols = lrows;
46345960306SStefano Zampini   } else { /* MatZeroRowsColumns implicitly assumes the rowscols indices are for a square matrix, here we handle a more general case */
46445960306SStefano Zampini     PetscInt i,nt,*t;
46545960306SStefano Zampini 
46645960306SStefano Zampini     ierr = PetscMalloc1(n,&t);CHKERRQ(ierr);
46745960306SStefano Zampini     for (i = 0, nt = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) t[nt++] = rowscols[i];
46845960306SStefano Zampini     ierr = PetscLayoutMapLocal(mat->cmap,nt,t,&nc,&lcols,NULL);CHKERRQ(ierr);
46945960306SStefano Zampini     ierr = PetscFree(t);CHKERRQ(ierr);
47045960306SStefano Zampini   }
47145960306SStefano Zampini   ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,nc,lcols,diag,PETSC_TRUE);CHKERRQ(ierr);
47245960306SStefano Zampini   if (!congruent) {
47345960306SStefano Zampini     ierr = PetscFree(lcols);CHKERRQ(ierr);
47445960306SStefano Zampini   }
47545960306SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
476ef5c7bd2SStefano Zampini   if (shell->axpy) {
477ef5c7bd2SStefano Zampini     ierr = MatZeroRowsColumns(shell->axpy,n,rowscols,0.0,NULL,NULL);CHKERRQ(ierr);
478ef5c7bd2SStefano Zampini   }
47945960306SStefano Zampini   PetscFunctionReturn(0);
48045960306SStefano Zampini }
48145960306SStefano Zampini 
482dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat)
483e51e0e81SBarry Smith {
484dfbe8321SBarry Smith   PetscErrorCode          ierr;
485bf0cc555SLisandro Dalcin   Mat_Shell               *shell = (Mat_Shell*)mat->data;
486b77ba244SStefano Zampini   MatShellMatFunctionList matmat;
487ed3cc1f0SBarry Smith 
4883a40ed3dSBarry Smith   PetscFunctionBegin;
48928f357e3SAlex Fikl   if (shell->ops->destroy) {
49028f357e3SAlex Fikl     ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr);
491bf0cc555SLisandro Dalcin   }
492ef5c7bd2SStefano Zampini   ierr = PetscMemzero(shell->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
4930c0fd78eSBarry Smith   ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
4940c0fd78eSBarry Smith   ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
4950c0fd78eSBarry Smith   ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
4968fe8eb27SJed Brown   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
4978fe8eb27SJed Brown   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
4985edf6516SJed Brown   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
4995edf6516SJed Brown   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
500b77ba244SStefano Zampini   ierr = VecDestroy(&shell->axpy_left);CHKERRQ(ierr);
501b77ba244SStefano Zampini   ierr = VecDestroy(&shell->axpy_right);CHKERRQ(ierr);
5029f137db4SBarry Smith   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
50345960306SStefano Zampini   ierr = VecDestroy(&shell->zvals_w);CHKERRQ(ierr);
50445960306SStefano Zampini   ierr = VecDestroy(&shell->zvals);CHKERRQ(ierr);
50545960306SStefano Zampini   ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr);
50645960306SStefano Zampini   ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr);
50745960306SStefano Zampini   ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr);
50845960306SStefano Zampini   ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr);
509b77ba244SStefano Zampini 
510b77ba244SStefano Zampini   matmat = shell->matmat;
511b77ba244SStefano Zampini   while (matmat) {
512b77ba244SStefano Zampini     MatShellMatFunctionList next = matmat->next;
513b77ba244SStefano Zampini 
514b77ba244SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)mat,matmat->composedname,NULL);CHKERRQ(ierr);
515b77ba244SStefano Zampini     ierr = PetscFree(matmat->composedname);CHKERRQ(ierr);
516b77ba244SStefano Zampini     ierr = PetscFree(matmat->resultname);CHKERRQ(ierr);
517b77ba244SStefano Zampini     ierr = PetscFree(matmat);CHKERRQ(ierr);
518b77ba244SStefano Zampini     matmat = next;
519b77ba244SStefano Zampini   }
520789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellGetContext_C",NULL);CHKERRQ(ierr);
521789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetContext_C",NULL);CHKERRQ(ierr);
522db77db73SJeremy L Thompson   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetVecType_C",NULL);CHKERRQ(ierr);
523789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetManageScalingShifts_C",NULL);CHKERRQ(ierr);
524789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetOperation_C",NULL);CHKERRQ(ierr);
525789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellGetOperation_C",NULL);CHKERRQ(ierr);
526b77ba244SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetMatProductOperation_C",NULL);CHKERRQ(ierr);
527b77ba244SStefano Zampini   ierr = PetscFree(mat->data);CHKERRQ(ierr);
528b77ba244SStefano Zampini   PetscFunctionReturn(0);
529b77ba244SStefano Zampini }
530b77ba244SStefano Zampini 
531b77ba244SStefano Zampini typedef struct {
532b77ba244SStefano Zampini   PetscErrorCode (*numeric)(Mat,Mat,Mat,void*);
533b77ba244SStefano Zampini   PetscErrorCode (*destroy)(void*);
534b77ba244SStefano Zampini   void           *userdata;
535b77ba244SStefano Zampini   Mat            B;
536b77ba244SStefano Zampini   Mat            Bt;
537b77ba244SStefano Zampini   Mat            axpy;
538b77ba244SStefano Zampini } MatMatDataShell;
539b77ba244SStefano Zampini 
540b77ba244SStefano Zampini static PetscErrorCode DestroyMatMatDataShell(void *data)
541b77ba244SStefano Zampini {
542b77ba244SStefano Zampini   MatMatDataShell *mmdata = (MatMatDataShell *)data;
543b77ba244SStefano Zampini   PetscErrorCode ierr;
544b77ba244SStefano Zampini 
545b77ba244SStefano Zampini   PetscFunctionBegin;
546b77ba244SStefano Zampini   if (mmdata->destroy) {
547b77ba244SStefano Zampini     ierr = (*mmdata->destroy)(mmdata->userdata);CHKERRQ(ierr);
548b77ba244SStefano Zampini   }
549b77ba244SStefano Zampini   ierr = MatDestroy(&mmdata->B);CHKERRQ(ierr);
550b77ba244SStefano Zampini   ierr = MatDestroy(&mmdata->Bt);CHKERRQ(ierr);
551b77ba244SStefano Zampini   ierr = MatDestroy(&mmdata->axpy);CHKERRQ(ierr);
552b77ba244SStefano Zampini   ierr = PetscFree(mmdata);CHKERRQ(ierr);
553b77ba244SStefano Zampini   PetscFunctionReturn(0);
554b77ba244SStefano Zampini }
555b77ba244SStefano Zampini 
556b77ba244SStefano Zampini static PetscErrorCode MatProductNumeric_Shell_X(Mat D)
557b77ba244SStefano Zampini {
558b77ba244SStefano Zampini   PetscErrorCode  ierr;
559b77ba244SStefano Zampini   Mat_Product     *product;
560b77ba244SStefano Zampini   Mat             A, B;
561b77ba244SStefano Zampini   MatMatDataShell *mdata;
562b77ba244SStefano Zampini   PetscScalar     zero = 0.0;
563b77ba244SStefano Zampini 
564b77ba244SStefano Zampini   PetscFunctionBegin;
565b77ba244SStefano Zampini   MatCheckProduct(D,1);
566b77ba244SStefano Zampini   product = D->product;
567b77ba244SStefano Zampini   if (!product->data) SETERRQ(PetscObjectComm((PetscObject)D),PETSC_ERR_PLIB,"Product data empty");
568b77ba244SStefano Zampini   A = product->A;
569b77ba244SStefano Zampini   B = product->B;
570b77ba244SStefano Zampini   mdata = (MatMatDataShell*)product->data;
571b77ba244SStefano Zampini   if (mdata->numeric) {
572b77ba244SStefano Zampini     Mat_Shell      *shell = (Mat_Shell*)A->data;
573b77ba244SStefano Zampini     PetscErrorCode (*stashsym)(Mat) = D->ops->productsymbolic;
574b77ba244SStefano Zampini     PetscErrorCode (*stashnum)(Mat) = D->ops->productnumeric;
575b77ba244SStefano Zampini     PetscBool      useBmdata = PETSC_FALSE, newB = PETSC_TRUE;
576b77ba244SStefano Zampini 
577b77ba244SStefano Zampini     if (shell->managescalingshifts) {
578b77ba244SStefano Zampini       if (shell->zcols || shell->zrows) SETERRQ(PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProduct not supported with zeroed rows/columns");
579b77ba244SStefano Zampini       if (shell->right || shell->left) {
580b77ba244SStefano Zampini         useBmdata = PETSC_TRUE;
581b77ba244SStefano Zampini         if (!mdata->B) {
582b77ba244SStefano Zampini           ierr = MatDuplicate(B,MAT_SHARE_NONZERO_PATTERN,&mdata->B);CHKERRQ(ierr);
583b77ba244SStefano Zampini         } else {
584b77ba244SStefano Zampini           newB = PETSC_FALSE;
585b77ba244SStefano Zampini         }
586b77ba244SStefano Zampini         ierr = MatCopy(B,mdata->B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
587b77ba244SStefano Zampini       }
588b77ba244SStefano Zampini       switch (product->type) {
589b77ba244SStefano Zampini       case MATPRODUCT_AB: /* s L A R B + v L R B + L D R B */
590b77ba244SStefano Zampini         if (shell->right) {
591b77ba244SStefano Zampini           ierr = MatDiagonalScale(mdata->B,shell->right,NULL);CHKERRQ(ierr);
592b77ba244SStefano Zampini         }
593b77ba244SStefano Zampini         break;
594b77ba244SStefano Zampini       case MATPRODUCT_AtB: /* s R A^t L B + v R L B + R D L B */
595b77ba244SStefano Zampini         if (shell->left) {
596b77ba244SStefano Zampini           ierr = MatDiagonalScale(mdata->B,shell->left,NULL);CHKERRQ(ierr);
597b77ba244SStefano Zampini         }
598b77ba244SStefano Zampini         break;
599b77ba244SStefano Zampini       case MATPRODUCT_ABt: /* s L A R B^t + v L R B^t + L D R B^t */
600b77ba244SStefano Zampini         if (shell->right) {
601b77ba244SStefano Zampini           ierr = MatDiagonalScale(mdata->B,NULL,shell->right);CHKERRQ(ierr);
602b77ba244SStefano Zampini         }
603b77ba244SStefano Zampini         break;
604b77ba244SStefano Zampini       case MATPRODUCT_RARt: /* s B L A R B^t + v B L R B^t + B L D R B^t */
605b77ba244SStefano Zampini         if (shell->right && shell->left) {
606b77ba244SStefano Zampini           PetscBool flg;
607b77ba244SStefano Zampini 
608b77ba244SStefano Zampini           ierr = VecEqual(shell->right,shell->left,&flg);CHKERRQ(ierr);
609b77ba244SStefano Zampini           if (!flg) SETERRQ3(PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProductSymbolic type %s not supported for %s and %s matrices because left scaling != from right scaling",MatProductTypes[product->type],((PetscObject)A)->type_name,((PetscObject)B)->type_name);
610b77ba244SStefano Zampini         }
611b77ba244SStefano Zampini         if (shell->right) {
612b77ba244SStefano Zampini           ierr = MatDiagonalScale(mdata->B,NULL,shell->right);CHKERRQ(ierr);
613b77ba244SStefano Zampini         }
614b77ba244SStefano Zampini         break;
615b77ba244SStefano Zampini       case MATPRODUCT_PtAP: /* s B^t L A R B + v B^t L R B + B^t L D R B */
616b77ba244SStefano Zampini         if (shell->right && shell->left) {
617b77ba244SStefano Zampini           PetscBool flg;
618b77ba244SStefano Zampini 
619b77ba244SStefano Zampini           ierr = VecEqual(shell->right,shell->left,&flg);CHKERRQ(ierr);
620b77ba244SStefano Zampini           if (!flg) SETERRQ3(PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProductSymbolic type %s not supported for %s and %s matrices because left scaling != from right scaling",MatProductTypes[product->type],((PetscObject)A)->type_name,((PetscObject)B)->type_name);
621b77ba244SStefano Zampini         }
622b77ba244SStefano Zampini         if (shell->right) {
623b77ba244SStefano Zampini           ierr = MatDiagonalScale(mdata->B,shell->right,NULL);CHKERRQ(ierr);
624b77ba244SStefano Zampini         }
625b77ba244SStefano Zampini         break;
626b77ba244SStefano Zampini       default: SETERRQ3(PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProductSymbolic type %s not supported for %s and %s matrices",MatProductTypes[product->type],((PetscObject)A)->type_name,((PetscObject)B)->type_name);
627b77ba244SStefano Zampini       }
628b77ba244SStefano Zampini     }
629b77ba244SStefano Zampini     /* allow the user to call MatMat operations on D */
630b77ba244SStefano Zampini     D->product = NULL;
631b77ba244SStefano Zampini     D->ops->productsymbolic = NULL;
632b77ba244SStefano Zampini     D->ops->productnumeric  = NULL;
633b77ba244SStefano Zampini 
634b77ba244SStefano Zampini     ierr = (*mdata->numeric)(A,useBmdata ? mdata->B : B,D,mdata->userdata);CHKERRQ(ierr);
635b77ba244SStefano Zampini 
636b77ba244SStefano Zampini     /* clear any leftover user data and restore D pointers */
637b77ba244SStefano Zampini     ierr = MatProductClear(D);CHKERRQ(ierr);
638b77ba244SStefano Zampini     D->ops->productsymbolic = stashsym;
639b77ba244SStefano Zampini     D->ops->productnumeric  = stashnum;
640b77ba244SStefano Zampini     D->product = product;
641b77ba244SStefano Zampini 
642b77ba244SStefano Zampini     if (shell->managescalingshifts) {
643b77ba244SStefano Zampini       ierr = MatScale(D,shell->vscale);CHKERRQ(ierr);
644b77ba244SStefano Zampini       switch (product->type) {
645b77ba244SStefano Zampini       case MATPRODUCT_AB: /* s L A R B + v L R B + L D R B */
646b77ba244SStefano Zampini       case MATPRODUCT_ABt: /* s L A R B^t + v L R B^t + L D R B^t */
647b77ba244SStefano Zampini         if (shell->left) {
648b77ba244SStefano Zampini           ierr = MatDiagonalScale(D,shell->left,NULL);CHKERRQ(ierr);
649b77ba244SStefano Zampini           if (shell->dshift || shell->vshift != zero) {
650b77ba244SStefano Zampini             if (!shell->left_work) {ierr = MatCreateVecs(A,NULL,&shell->left_work);CHKERRQ(ierr);}
651b77ba244SStefano Zampini             if (shell->dshift) {
652b77ba244SStefano Zampini               ierr = VecCopy(shell->dshift,shell->left_work);CHKERRQ(ierr);
653b77ba244SStefano Zampini               ierr = VecShift(shell->left_work,shell->vshift);CHKERRQ(ierr);
654b77ba244SStefano Zampini               ierr = VecPointwiseMult(shell->left_work,shell->left_work,shell->left);CHKERRQ(ierr);
655b77ba244SStefano Zampini             } else {
656b77ba244SStefano Zampini               ierr = VecSet(shell->left_work,shell->vshift);CHKERRQ(ierr);
657b77ba244SStefano Zampini             }
658b77ba244SStefano Zampini             if (product->type == MATPRODUCT_ABt) {
659b77ba244SStefano Zampini               MatReuse     reuse = mdata->Bt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
660b77ba244SStefano Zampini               MatStructure str = mdata->Bt ? SUBSET_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
661b77ba244SStefano Zampini 
662b77ba244SStefano Zampini               ierr = MatTranspose(mdata->B,reuse,&mdata->Bt);CHKERRQ(ierr);
663b77ba244SStefano Zampini               ierr = MatDiagonalScale(mdata->Bt,shell->left_work,NULL);CHKERRQ(ierr);
664b77ba244SStefano Zampini               ierr = MatAXPY(D,1.0,mdata->Bt,str);CHKERRQ(ierr);
665b77ba244SStefano Zampini             } else {
666b77ba244SStefano Zampini               MatStructure str = newB ? DIFFERENT_NONZERO_PATTERN : SUBSET_NONZERO_PATTERN;
667b77ba244SStefano Zampini 
668b77ba244SStefano Zampini               ierr = MatDiagonalScale(mdata->B,shell->left_work,NULL);CHKERRQ(ierr);
669b77ba244SStefano Zampini               ierr = MatAXPY(D,1.0,mdata->B,str);CHKERRQ(ierr);
670b77ba244SStefano Zampini             }
671b77ba244SStefano Zampini           }
672b77ba244SStefano Zampini         }
673b77ba244SStefano Zampini         break;
674b77ba244SStefano Zampini       case MATPRODUCT_AtB: /* s R A^t L B + v R L B + R D L B */
675b77ba244SStefano Zampini         if (shell->right) {
676b77ba244SStefano Zampini           ierr = MatDiagonalScale(D,shell->right,NULL);CHKERRQ(ierr);
677b77ba244SStefano Zampini           if (shell->dshift || shell->vshift != zero) {
678b77ba244SStefano Zampini             MatStructure str = newB ? DIFFERENT_NONZERO_PATTERN : SUBSET_NONZERO_PATTERN;
679b77ba244SStefano Zampini 
680b77ba244SStefano Zampini             if (!shell->right_work) {ierr = MatCreateVecs(A,&shell->right_work,NULL);CHKERRQ(ierr);}
681b77ba244SStefano Zampini             if (shell->dshift) {
682b77ba244SStefano Zampini               ierr = VecCopy(shell->dshift,shell->right_work);CHKERRQ(ierr);
683b77ba244SStefano Zampini               ierr = VecShift(shell->right_work,shell->vshift);CHKERRQ(ierr);
684b77ba244SStefano Zampini               ierr = VecPointwiseMult(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
685b77ba244SStefano Zampini             } else {
686b77ba244SStefano Zampini               ierr = VecSet(shell->right_work,shell->vshift);CHKERRQ(ierr);
687b77ba244SStefano Zampini             }
688b77ba244SStefano Zampini             ierr = MatDiagonalScale(mdata->B,shell->right_work,NULL);CHKERRQ(ierr);
689b77ba244SStefano Zampini             ierr = MatAXPY(D,1.0,mdata->B,str);CHKERRQ(ierr);
690b77ba244SStefano Zampini           }
691b77ba244SStefano Zampini         }
692b77ba244SStefano Zampini         break;
693b77ba244SStefano Zampini       case MATPRODUCT_PtAP: /* s B^t L A R B + v B^t L R B + B^t L D R B */
694b77ba244SStefano Zampini       case MATPRODUCT_RARt: /* s B L A R B^t + v B L R B^t + B L D R B^t */
695b77ba244SStefano Zampini         if (shell->dshift || shell->vshift != zero) SETERRQ3(PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProductSymbolic type %s not supported for %s and %s matrices with diagonal shift",MatProductTypes[product->type],((PetscObject)A)->type_name,((PetscObject)B)->type_name);
696b77ba244SStefano Zampini         break;
697b77ba244SStefano Zampini       default: SETERRQ3(PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProductSymbolic type %s not supported for %s and %s matrices",MatProductTypes[product->type],((PetscObject)A)->type_name,((PetscObject)B)->type_name);
698b77ba244SStefano Zampini       }
699b77ba244SStefano Zampini       if (shell->axpy && shell->axpy_vscale != zero) {
700b77ba244SStefano Zampini         Mat              X;
701b77ba244SStefano Zampini         PetscObjectState axpy_state;
702b77ba244SStefano Zampini         MatStructure     str = DIFFERENT_NONZERO_PATTERN; /* not sure it is safe to ever use SUBSET_NONZERO_PATTERN */
703b77ba244SStefano Zampini 
704b77ba244SStefano Zampini         ierr = MatShellGetContext(shell->axpy,(void *)&X);CHKERRQ(ierr);
705b77ba244SStefano Zampini         ierr = PetscObjectStateGet((PetscObject)X,&axpy_state);CHKERRQ(ierr);
706b77ba244SStefano Zampini         if (shell->axpy_state != axpy_state) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)");
707b77ba244SStefano Zampini         if (!mdata->axpy) {
708b77ba244SStefano Zampini           str  = DIFFERENT_NONZERO_PATTERN;
709b77ba244SStefano Zampini           ierr = MatProductCreate(shell->axpy,B,NULL,&mdata->axpy);CHKERRQ(ierr);
710b77ba244SStefano Zampini           ierr = MatProductSetType(mdata->axpy,product->type);CHKERRQ(ierr);
711b77ba244SStefano Zampini           ierr = MatProductSetFromOptions(mdata->axpy);CHKERRQ(ierr);
712b77ba244SStefano Zampini           ierr = MatProductSymbolic(mdata->axpy);CHKERRQ(ierr);
713b77ba244SStefano Zampini         } else { /* May be that shell->axpy has changed */
714b77ba244SStefano Zampini           PetscBool flg;
715b77ba244SStefano Zampini 
716b77ba244SStefano Zampini           ierr = MatProductReplaceMats(shell->axpy,B,NULL,mdata->axpy);CHKERRQ(ierr);
717b77ba244SStefano Zampini           ierr = MatHasOperation(mdata->axpy,MATOP_PRODUCTSYMBOLIC,&flg);CHKERRQ(ierr);
718b77ba244SStefano Zampini           if (!flg) {
719b77ba244SStefano Zampini             str  = DIFFERENT_NONZERO_PATTERN;
720b77ba244SStefano Zampini             ierr = MatProductSetFromOptions(mdata->axpy);CHKERRQ(ierr);
721b77ba244SStefano Zampini             ierr = MatProductSymbolic(mdata->axpy);CHKERRQ(ierr);
722b77ba244SStefano Zampini           }
723b77ba244SStefano Zampini         }
724b77ba244SStefano Zampini         ierr = MatProductNumeric(mdata->axpy);CHKERRQ(ierr);
725b77ba244SStefano Zampini         ierr = MatAXPY(D,shell->axpy_vscale,mdata->axpy,str);CHKERRQ(ierr);
726b77ba244SStefano Zampini       }
727b77ba244SStefano Zampini     }
728b77ba244SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)D),PETSC_ERR_PLIB,"Missing numeric operation");
729b77ba244SStefano Zampini   PetscFunctionReturn(0);
730b77ba244SStefano Zampini }
731b77ba244SStefano Zampini 
732b77ba244SStefano Zampini static PetscErrorCode MatProductSymbolic_Shell_X(Mat D)
733b77ba244SStefano Zampini {
734b77ba244SStefano Zampini   PetscErrorCode          ierr;
735b77ba244SStefano Zampini   Mat_Product             *product;
736b77ba244SStefano Zampini   Mat                     A,B;
737b77ba244SStefano Zampini   MatShellMatFunctionList matmat;
738b77ba244SStefano Zampini   Mat_Shell               *shell;
739b77ba244SStefano Zampini   PetscBool               flg;
740b77ba244SStefano Zampini   char                    composedname[256];
741b77ba244SStefano Zampini   MatMatDataShell         *mdata;
742b77ba244SStefano Zampini 
743b77ba244SStefano Zampini   PetscFunctionBegin;
744b77ba244SStefano Zampini   MatCheckProduct(D,1);
745b77ba244SStefano Zampini   product = D->product;
746b77ba244SStefano Zampini   if (product->data) SETERRQ(PetscObjectComm((PetscObject)D),PETSC_ERR_PLIB,"Product data not empty");
747b77ba244SStefano Zampini   A = product->A;
748b77ba244SStefano Zampini   B = product->B;
749b77ba244SStefano Zampini   shell = (Mat_Shell*)A->data;
750b77ba244SStefano Zampini   matmat = shell->matmat;
751b77ba244SStefano Zampini   ierr = PetscSNPrintf(composedname,sizeof(composedname),"MatProductSetFromOptions_%s_%s_C",((PetscObject)A)->type_name,((PetscObject)B)->type_name);CHKERRQ(ierr);
752b77ba244SStefano Zampini   while (matmat) {
753b77ba244SStefano Zampini     ierr = PetscStrcmp(composedname,matmat->composedname,&flg);CHKERRQ(ierr);
754b77ba244SStefano Zampini     flg  = (PetscBool)(flg && (matmat->ptype == product->type));
755b77ba244SStefano Zampini     if (flg) break;
756b77ba244SStefano Zampini     matmat = matmat->next;
757b77ba244SStefano Zampini   }
758b77ba244SStefano Zampini   if (!flg) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Composedname \"%s\" for product type %s not found",composedname,MatProductTypes[product->type]);
759b77ba244SStefano Zampini   switch (product->type) {
760b77ba244SStefano Zampini   case MATPRODUCT_AB:
761b77ba244SStefano Zampini     ierr = MatSetSizes(D,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
762b77ba244SStefano Zampini     break;
763b77ba244SStefano Zampini   case MATPRODUCT_AtB:
764b77ba244SStefano Zampini     ierr = MatSetSizes(D,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N);CHKERRQ(ierr);
765b77ba244SStefano Zampini     break;
766b77ba244SStefano Zampini   case MATPRODUCT_ABt:
767b77ba244SStefano Zampini     ierr = MatSetSizes(D,A->rmap->n,B->rmap->n,A->rmap->N,B->rmap->N);CHKERRQ(ierr);
768b77ba244SStefano Zampini     break;
769b77ba244SStefano Zampini   case MATPRODUCT_RARt:
770b77ba244SStefano Zampini     ierr = MatSetSizes(D,B->rmap->n,B->rmap->n,B->rmap->N,B->rmap->N);CHKERRQ(ierr);
771b77ba244SStefano Zampini     break;
772b77ba244SStefano Zampini   case MATPRODUCT_PtAP:
773b77ba244SStefano Zampini     ierr = MatSetSizes(D,B->cmap->n,B->cmap->n,B->cmap->N,B->cmap->N);CHKERRQ(ierr);
774b77ba244SStefano Zampini     break;
775b77ba244SStefano Zampini   default: SETERRQ3(PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProductSymbolic type %s not supported for %s and %s matrices",MatProductTypes[product->type],((PetscObject)A)->type_name,((PetscObject)B)->type_name);
776b77ba244SStefano Zampini   }
777b77ba244SStefano Zampini   /* respect users who passed in a matrix for which resultname is the base type */
778b77ba244SStefano Zampini   if (matmat->resultname) {
779b77ba244SStefano Zampini     ierr = PetscObjectBaseTypeCompare((PetscObject)D,matmat->resultname,&flg);CHKERRQ(ierr);
780b77ba244SStefano Zampini     if (!flg) {
781b77ba244SStefano Zampini       ierr = MatSetType(D,matmat->resultname);CHKERRQ(ierr);
782b77ba244SStefano Zampini     }
783b77ba244SStefano Zampini   }
784b77ba244SStefano Zampini   /* If matrix type was not set or different, we need to reset this pointers */
785b77ba244SStefano Zampini   D->ops->productsymbolic = MatProductSymbolic_Shell_X;
786b77ba244SStefano Zampini   D->ops->productnumeric  = MatProductNumeric_Shell_X;
787b77ba244SStefano Zampini   /* attach product data */
788b77ba244SStefano Zampini   ierr = PetscNew(&mdata);CHKERRQ(ierr);
789b77ba244SStefano Zampini   mdata->numeric = matmat->numeric;
790b77ba244SStefano Zampini   mdata->destroy = matmat->destroy;
791b77ba244SStefano Zampini   if (matmat->symbolic) {
792b77ba244SStefano Zampini     ierr = (*matmat->symbolic)(A,B,D,&mdata->userdata);CHKERRQ(ierr);
793b77ba244SStefano Zampini   } else { /* call general setup if symbolic operation not provided */
794b77ba244SStefano Zampini     ierr = MatSetUp(D);CHKERRQ(ierr);
795b77ba244SStefano Zampini   }
796b77ba244SStefano Zampini   if (!D->product) SETERRQ(PetscObjectComm((PetscObject)D),PETSC_ERR_COR,"Product disappeared after user symbolic phase");
797b77ba244SStefano Zampini   if (D->product->data) SETERRQ(PetscObjectComm((PetscObject)D),PETSC_ERR_COR,"Product data not empty after user symbolic phase");
798b77ba244SStefano Zampini   D->product->data = mdata;
799b77ba244SStefano Zampini   D->product->destroy = DestroyMatMatDataShell;
800b77ba244SStefano Zampini   /* Be sure to reset these pointers if the user did something unexpected */
801b77ba244SStefano Zampini   D->ops->productsymbolic = MatProductSymbolic_Shell_X;
802b77ba244SStefano Zampini   D->ops->productnumeric  = MatProductNumeric_Shell_X;
803b77ba244SStefano Zampini   PetscFunctionReturn(0);
804b77ba244SStefano Zampini }
805b77ba244SStefano Zampini 
806b77ba244SStefano Zampini static PetscErrorCode MatProductSetFromOptions_Shell_X(Mat D)
807b77ba244SStefano Zampini {
808b77ba244SStefano Zampini   PetscErrorCode          ierr;
809b77ba244SStefano Zampini   Mat_Product             *product;
810b77ba244SStefano Zampini   Mat                     A,B;
811b77ba244SStefano Zampini   MatShellMatFunctionList matmat;
812b77ba244SStefano Zampini   Mat_Shell               *shell;
813b77ba244SStefano Zampini   PetscBool               flg;
814b77ba244SStefano Zampini   char                    composedname[256];
815b77ba244SStefano Zampini 
816b77ba244SStefano Zampini   PetscFunctionBegin;
817b77ba244SStefano Zampini   MatCheckProduct(D,1);
818b77ba244SStefano Zampini   product = D->product;
819b77ba244SStefano Zampini   A = product->A;
820b77ba244SStefano Zampini   B = product->B;
821b77ba244SStefano Zampini   ierr = MatIsShell(A,&flg);CHKERRQ(ierr);
822b77ba244SStefano Zampini   if (!flg) PetscFunctionReturn(0);
823b77ba244SStefano Zampini   shell = (Mat_Shell*)A->data;
824b77ba244SStefano Zampini   matmat = shell->matmat;
825b77ba244SStefano Zampini   ierr = PetscSNPrintf(composedname,sizeof(composedname),"MatProductSetFromOptions_%s_%s_C",((PetscObject)A)->type_name,((PetscObject)B)->type_name);CHKERRQ(ierr);
826b77ba244SStefano Zampini   while (matmat) {
827b77ba244SStefano Zampini     ierr = PetscStrcmp(composedname,matmat->composedname,&flg);CHKERRQ(ierr);
828b77ba244SStefano Zampini     flg  = (PetscBool)(flg && (matmat->ptype == product->type));
829b77ba244SStefano Zampini     if (flg) break;
830b77ba244SStefano Zampini     matmat = matmat->next;
831b77ba244SStefano Zampini   }
832b77ba244SStefano Zampini   if (flg) { D->ops->productsymbolic = MatProductSymbolic_Shell_X; }
833b77ba244SStefano Zampini   else { ierr = PetscInfo2(D,"  symbolic product %s not registered for product type %s\n",composedname,MatProductTypes[product->type]);CHKERRQ(ierr); }
834b77ba244SStefano Zampini   PetscFunctionReturn(0);
835b77ba244SStefano Zampini }
836b77ba244SStefano Zampini 
837b77ba244SStefano Zampini static PetscErrorCode MatShellSetMatProductOperation_Private(Mat A,MatProductType ptype,PetscErrorCode (*symbolic)(Mat,Mat,Mat,void**),PetscErrorCode (*numeric)(Mat,Mat,Mat,void*),PetscErrorCode (*destroy)(void*),char *composedname,const char *resultname)
838b77ba244SStefano Zampini {
839b77ba244SStefano Zampini   PetscBool               flg;
840b77ba244SStefano Zampini   PetscErrorCode          ierr;
841b77ba244SStefano Zampini   Mat_Shell               *shell;
842b77ba244SStefano Zampini   MatShellMatFunctionList matmat;
843b77ba244SStefano Zampini 
844b77ba244SStefano Zampini   PetscFunctionBegin;
845b77ba244SStefano Zampini   if (!numeric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing numeric routine");
846b77ba244SStefano Zampini   if (!composedname) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing composed name");
847b77ba244SStefano Zampini 
848b77ba244SStefano Zampini   /* add product callback */
849b77ba244SStefano Zampini   shell = (Mat_Shell*)A->data;
850b77ba244SStefano Zampini   matmat = shell->matmat;
851b77ba244SStefano Zampini   if (!matmat) {
852b77ba244SStefano Zampini     ierr = PetscNew(&shell->matmat);CHKERRQ(ierr);
853b77ba244SStefano Zampini     matmat = shell->matmat;
854b77ba244SStefano Zampini   } else {
855b77ba244SStefano Zampini     MatShellMatFunctionList entry = matmat;
856b77ba244SStefano Zampini     while (entry) {
857b77ba244SStefano Zampini       ierr = PetscStrcmp(composedname,entry->composedname,&flg);CHKERRQ(ierr);
858b77ba244SStefano Zampini       flg  = (PetscBool)(flg && (entry->ptype == ptype));
859*843d480fSStefano Zampini       if (flg) goto set;
860b77ba244SStefano Zampini       matmat = entry;
861b77ba244SStefano Zampini       entry = entry->next;
862b77ba244SStefano Zampini     }
863b77ba244SStefano Zampini     ierr = PetscNew(&matmat->next);CHKERRQ(ierr);
864b77ba244SStefano Zampini     matmat = matmat->next;
865b77ba244SStefano Zampini   }
866b77ba244SStefano Zampini 
867*843d480fSStefano Zampini set:
868b77ba244SStefano Zampini   matmat->symbolic = symbolic;
869b77ba244SStefano Zampini   matmat->numeric  = numeric;
870b77ba244SStefano Zampini   matmat->destroy  = destroy;
871b77ba244SStefano Zampini   matmat->ptype    = ptype;
872b77ba244SStefano Zampini   ierr = PetscFree(matmat->composedname);CHKERRQ(ierr);
873b77ba244SStefano Zampini   ierr = PetscFree(matmat->resultname);CHKERRQ(ierr);
874b77ba244SStefano Zampini   ierr = PetscStrallocpy(composedname,&matmat->composedname);CHKERRQ(ierr);
875b77ba244SStefano Zampini   ierr = PetscStrallocpy(resultname,&matmat->resultname);CHKERRQ(ierr);
876b77ba244SStefano Zampini   ierr = PetscInfo3(A,"Composing %s for product type %s with result %s\n",matmat->composedname,MatProductTypes[matmat->ptype],matmat->resultname ? matmat->resultname : "not specified");CHKERRQ(ierr);
877b77ba244SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,matmat->composedname,MatProductSetFromOptions_Shell_X);CHKERRQ(ierr);
878b77ba244SStefano Zampini   PetscFunctionReturn(0);
879b77ba244SStefano Zampini }
880b77ba244SStefano Zampini 
881b77ba244SStefano Zampini /*@C
882b77ba244SStefano Zampini     MatShellSetMatProductOperation - Allows user to set a matrix matrix operation for a shell matrix.
883b77ba244SStefano Zampini 
884b77ba244SStefano Zampini    Logically Collective on Mat
885b77ba244SStefano Zampini 
886b77ba244SStefano Zampini     Input Parameters:
887b77ba244SStefano Zampini +   A - the shell matrix
888b77ba244SStefano Zampini .   ptype - the product type
889b77ba244SStefano Zampini .   symbolic - the function for the symbolic phase (can be NULL)
890b77ba244SStefano Zampini .   numeric - the function for the numerical phase
891b77ba244SStefano Zampini .   destroy - the function for the destruction of the needed data generated during the symbolic phase (can be NULL)
892b77ba244SStefano Zampini .   Btype - the matrix type for the matrix to be multiplied against
893b77ba244SStefano Zampini -   Ctype - the matrix type for the result (can be NULL)
894b77ba244SStefano Zampini 
895b77ba244SStefano Zampini    Level: advanced
896b77ba244SStefano Zampini 
897b77ba244SStefano Zampini     Usage:
898b77ba244SStefano Zampini $      extern PetscErrorCode usersymbolic(Mat,Mat,Mat,void**);
899b77ba244SStefano Zampini $      extern PetscErrorCode usernumeric(Mat,Mat,Mat,void*);
900b77ba244SStefano Zampini $      extern PetscErrorCode userdestroy(void*);
901b77ba244SStefano Zampini $      MatCreateShell(comm,m,n,M,N,ctx,&A);
902b77ba244SStefano Zampini $      MatShellSetMatProductOperation(A,MATPRODUCT_AB,usersymbolic,usernumeric,userdestroy,MATSEQAIJ,MATDENSE);
903b77ba244SStefano Zampini $      [ create B of type SEQAIJ etc..]
904b77ba244SStefano Zampini $      MatProductCreate(A,B,NULL,&C);
905b77ba244SStefano Zampini $      MatProductSetType(C,MATPRODUCT_AB);
906b77ba244SStefano Zampini $      MatProductSetFromOptions(C);
907b77ba244SStefano Zampini $      MatProductSymbolic(C); -> actually runs the user defined symbolic operation
908b77ba244SStefano Zampini $      MatProductNumeric(C); -> actually runs the user defined numeric operation
909b77ba244SStefano Zampini $      [ use C = A*B ]
910b77ba244SStefano Zampini 
911b77ba244SStefano Zampini     Notes:
912b77ba244SStefano Zampini     MATPRODUCT_ABC is not supported yet. Not supported in Fortran.
913b77ba244SStefano Zampini     If the symbolic phase is not specified, MatSetUp() is called on the result matrix that must have its type set if Ctype is NULL.
914b77ba244SStefano Zampini     Any additional data needed by the matrix product needs to be returned during the symbolic phase and destroyed with the destroy callback.
915b77ba244SStefano Zampini     PETSc will take care of calling the user-defined callbacks.
916b77ba244SStefano Zampini     It is allowed to specify the same callbacks for different Btype matrix types.
917b77ba244SStefano Zampini     The couple (Btype,ptype) uniquely identifies the operation: the last specified callbacks takes precedence.
918b77ba244SStefano Zampini 
919b77ba244SStefano Zampini .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatProductType, MatType, MatSetUp()
920b77ba244SStefano Zampini @*/
921b77ba244SStefano Zampini PetscErrorCode MatShellSetMatProductOperation(Mat A,MatProductType ptype,PetscErrorCode (*symbolic)(Mat,Mat,Mat,void**),PetscErrorCode (*numeric)(Mat,Mat,Mat,void*),PetscErrorCode (*destroy)(void *),MatType Btype,MatType Ctype)
922b77ba244SStefano Zampini {
923b77ba244SStefano Zampini   PetscErrorCode ierr;
924b77ba244SStefano Zampini 
925b77ba244SStefano Zampini   PetscFunctionBegin;
926b77ba244SStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
927b77ba244SStefano Zampini   PetscValidLogicalCollectiveEnum(A,ptype,2);
928b77ba244SStefano Zampini   if (ptype == MATPRODUCT_ABC) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not for product type %s",MatProductTypes[ptype]);
929b77ba244SStefano Zampini   if (!numeric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Missing numeric routine, argument 4");
930b77ba244SStefano Zampini   PetscValidPointer(Btype,6);
931b77ba244SStefano Zampini   if (Ctype) PetscValidPointer(Ctype,7);
932b77ba244SStefano Zampini   ierr = PetscTryMethod(A,"MatShellSetMatProductOperation_C",(Mat,MatProductType,PetscErrorCode(*)(Mat,Mat,Mat,void**),PetscErrorCode(*)(Mat,Mat,Mat,void*),PetscErrorCode(*)(void*),MatType,MatType),(A,ptype,symbolic,numeric,destroy,Btype,Ctype));CHKERRQ(ierr);
933b77ba244SStefano Zampini   PetscFunctionReturn(0);
934b77ba244SStefano Zampini }
935b77ba244SStefano Zampini 
936b77ba244SStefano Zampini PetscErrorCode MatShellSetMatProductOperation_Shell(Mat A,MatProductType ptype,PetscErrorCode (*symbolic)(Mat,Mat,Mat,void**),PetscErrorCode (*numeric)(Mat,Mat,Mat,void*),PetscErrorCode (*destroy)(void *),MatType Btype,MatType Ctype)
937b77ba244SStefano Zampini {
938b77ba244SStefano Zampini   PetscBool      flg;
939b77ba244SStefano Zampini   PetscErrorCode ierr;
940b77ba244SStefano Zampini   char           composedname[256];
941b77ba244SStefano Zampini   MatRootName    Bnames = MatRootNameList, Cnames = MatRootNameList;
942b77ba244SStefano Zampini   PetscMPIInt    size;
943b77ba244SStefano Zampini 
944b77ba244SStefano Zampini   PetscFunctionBegin;
945b77ba244SStefano Zampini   PetscValidType(A,1);
946b77ba244SStefano Zampini   while (Bnames) { /* user passed in the root name */
947b77ba244SStefano Zampini     ierr = PetscStrcmp(Btype,Bnames->rname,&flg);CHKERRQ(ierr);
948b77ba244SStefano Zampini     if (flg) break;
949b77ba244SStefano Zampini     Bnames = Bnames->next;
950b77ba244SStefano Zampini   }
951b77ba244SStefano Zampini   while (Cnames) { /* user passed in the root name */
952b77ba244SStefano Zampini     ierr = PetscStrcmp(Ctype,Cnames->rname,&flg);CHKERRQ(ierr);
953b77ba244SStefano Zampini     if (flg) break;
954b77ba244SStefano Zampini     Cnames = Cnames->next;
955b77ba244SStefano Zampini   }
956ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr);
957b77ba244SStefano Zampini   Btype = Bnames ? (size > 1 ? Bnames->mname : Bnames->sname) : Btype;
958b77ba244SStefano Zampini   Ctype = Cnames ? (size > 1 ? Cnames->mname : Cnames->sname) : Ctype;
959b77ba244SStefano Zampini   ierr = PetscSNPrintf(composedname,sizeof(composedname),"MatProductSetFromOptions_%s_%s_C",((PetscObject)A)->type_name,Btype);CHKERRQ(ierr);
960b77ba244SStefano Zampini   ierr = MatShellSetMatProductOperation_Private(A,ptype,symbolic,numeric,destroy,composedname,Ctype);CHKERRQ(ierr);
9613a40ed3dSBarry Smith   PetscFunctionReturn(0);
962e51e0e81SBarry Smith }
963e51e0e81SBarry Smith 
9647fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
9657fabda3fSAlex Fikl {
96628f357e3SAlex Fikl   Mat_Shell               *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
9677fabda3fSAlex Fikl   PetscErrorCode          ierr;
968cb8c8a77SBarry Smith   PetscBool               matflg;
969b77ba244SStefano Zampini   MatShellMatFunctionList matmatA;
9707fabda3fSAlex Fikl 
9717fabda3fSAlex Fikl   PetscFunctionBegin;
972b77ba244SStefano Zampini   ierr = MatIsShell(B,&matflg);CHKERRQ(ierr);
973b77ba244SStefano Zampini   if (!matflg) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix %s not derived from MATSHELL",((PetscObject)B)->type_name);
9747fabda3fSAlex Fikl 
975cb8c8a77SBarry Smith   ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
976cb8c8a77SBarry Smith   ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
9777fabda3fSAlex Fikl 
978cb8c8a77SBarry Smith   if (shellA->ops->copy) {
97928f357e3SAlex Fikl     ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr);
980cb8c8a77SBarry Smith   }
9817fabda3fSAlex Fikl   shellB->vscale = shellA->vscale;
9827fabda3fSAlex Fikl   shellB->vshift = shellA->vshift;
9830c0fd78eSBarry Smith   if (shellA->dshift) {
9840c0fd78eSBarry Smith     if (!shellB->dshift) {
9850c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr);
9867fabda3fSAlex Fikl     }
9870c0fd78eSBarry Smith     ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr);
9887fabda3fSAlex Fikl   } else {
9899f137db4SBarry Smith     ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr);
9907fabda3fSAlex Fikl   }
9910c0fd78eSBarry Smith   if (shellA->left) {
9920c0fd78eSBarry Smith     if (!shellB->left) {
9930c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr);
9947fabda3fSAlex Fikl     }
9950c0fd78eSBarry Smith     ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr);
9967fabda3fSAlex Fikl   } else {
9979f137db4SBarry Smith     ierr = VecDestroy(&shellB->left);CHKERRQ(ierr);
9987fabda3fSAlex Fikl   }
9990c0fd78eSBarry Smith   if (shellA->right) {
10000c0fd78eSBarry Smith     if (!shellB->right) {
10010c0fd78eSBarry Smith       ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr);
10027fabda3fSAlex Fikl     }
10030c0fd78eSBarry Smith     ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr);
10047fabda3fSAlex Fikl   } else {
10059f137db4SBarry Smith     ierr = VecDestroy(&shellB->right);CHKERRQ(ierr);
10067fabda3fSAlex Fikl   }
100740e381c3SBarry Smith   ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr);
1008ef5c7bd2SStefano Zampini   shellB->axpy_vscale = 0.0;
1009ef5c7bd2SStefano Zampini   shellB->axpy_state  = 0;
101040e381c3SBarry Smith   if (shellA->axpy) {
101140e381c3SBarry Smith     ierr                = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr);
101240e381c3SBarry Smith     shellB->axpy        = shellA->axpy;
101340e381c3SBarry Smith     shellB->axpy_vscale = shellA->axpy_vscale;
1014ef5c7bd2SStefano Zampini     shellB->axpy_state  = shellA->axpy_state;
101540e381c3SBarry Smith   }
101645960306SStefano Zampini   if (shellA->zrows) {
101745960306SStefano Zampini     ierr = ISDuplicate(shellA->zrows,&shellB->zrows);CHKERRQ(ierr);
101845960306SStefano Zampini     if (shellA->zcols) {
101945960306SStefano Zampini       ierr = ISDuplicate(shellA->zcols,&shellB->zcols);CHKERRQ(ierr);
102045960306SStefano Zampini     }
102145960306SStefano Zampini     ierr = VecDuplicate(shellA->zvals,&shellB->zvals);CHKERRQ(ierr);
102245960306SStefano Zampini     ierr = VecCopy(shellA->zvals,shellB->zvals);CHKERRQ(ierr);
102345960306SStefano Zampini     ierr = VecDuplicate(shellA->zvals_w,&shellB->zvals_w);CHKERRQ(ierr);
102445960306SStefano Zampini     ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_r);CHKERRQ(ierr);
102545960306SStefano Zampini     ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_c);CHKERRQ(ierr);
102645960306SStefano Zampini     shellB->zvals_sct_r = shellA->zvals_sct_r;
102745960306SStefano Zampini     shellB->zvals_sct_c = shellA->zvals_sct_c;
102845960306SStefano Zampini   }
1029b77ba244SStefano Zampini 
1030b77ba244SStefano Zampini   matmatA = shellA->matmat;
1031b77ba244SStefano Zampini   if (matmatA) {
1032b77ba244SStefano Zampini     while (matmatA->next) {
1033b77ba244SStefano Zampini       ierr = MatShellSetMatProductOperation_Private(B,matmatA->ptype,matmatA->symbolic,matmatA->numeric,matmatA->destroy,matmatA->composedname,matmatA->resultname);CHKERRQ(ierr);
1034b77ba244SStefano Zampini       matmatA = matmatA->next;
1035b77ba244SStefano Zampini     }
1036b77ba244SStefano Zampini   }
10377fabda3fSAlex Fikl   PetscFunctionReturn(0);
10387fabda3fSAlex Fikl }
10397fabda3fSAlex Fikl 
1040cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M)
1041cb8c8a77SBarry Smith {
1042cb8c8a77SBarry Smith   PetscErrorCode ierr;
1043cb8c8a77SBarry Smith   void           *ctx;
1044cb8c8a77SBarry Smith 
1045cb8c8a77SBarry Smith   PetscFunctionBegin;
1046cb8c8a77SBarry Smith   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
1047cb8c8a77SBarry Smith   ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr);
1048b77ba244SStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)(*M),((PetscObject)mat)->type_name);CHKERRQ(ierr);
1049a4b1380bSStefano Zampini   if (op != MAT_DO_NOT_COPY_VALUES) {
1050cb8c8a77SBarry Smith     ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1051a4b1380bSStefano Zampini   }
1052cb8c8a77SBarry Smith   PetscFunctionReturn(0);
1053cb8c8a77SBarry Smith }
1054cb8c8a77SBarry Smith 
1055dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
1056ef66eb69SBarry Smith {
1057ef66eb69SBarry Smith   Mat_Shell        *shell = (Mat_Shell*)A->data;
1058dfbe8321SBarry Smith   PetscErrorCode   ierr;
105925578ef6SJed Brown   Vec              xx;
1060e3f487b0SBarry Smith   PetscObjectState instate,outstate;
1061ef66eb69SBarry Smith 
1062ef66eb69SBarry Smith   PetscFunctionBegin;
1063976e8c5aSLisandro Dalcin   if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL");
106445960306SStefano Zampini   ierr = MatShellPreZeroRight(A,x,&xx);CHKERRQ(ierr);
106545960306SStefano Zampini   ierr = MatShellPreScaleRight(A,xx,&xx);CHKERRQ(ierr);
106645960306SStefano Zampini   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
106728f357e3SAlex Fikl   ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr);
1068e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
1069e3f487b0SBarry Smith   if (instate == outstate) {
1070e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
1071e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
1072e3f487b0SBarry Smith   }
10738fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
10748fe8eb27SJed Brown   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
107545960306SStefano Zampini   ierr = MatShellPostZeroLeft(A,y);CHKERRQ(ierr);
10769f137db4SBarry Smith 
10779f137db4SBarry Smith   if (shell->axpy) {
1078ef5c7bd2SStefano Zampini     Mat              X;
1079ef5c7bd2SStefano Zampini     PetscObjectState axpy_state;
1080ef5c7bd2SStefano Zampini 
1081ef5c7bd2SStefano Zampini     ierr = MatShellGetContext(shell->axpy,(void *)&X);CHKERRQ(ierr);
1082ef5c7bd2SStefano Zampini     ierr = PetscObjectStateGet((PetscObject)X,&axpy_state);CHKERRQ(ierr);
108331c70108SStefano Zampini     if (shell->axpy_state != axpy_state) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)");
1084b77ba244SStefano Zampini 
1085b77ba244SStefano Zampini     ierr = MatCreateVecs(shell->axpy,shell->axpy_right ? NULL : &shell->axpy_right,shell->axpy_left ? NULL : &shell->axpy_left);CHKERRQ(ierr);
1086b77ba244SStefano Zampini     ierr = VecCopy(x,shell->axpy_right);CHKERRQ(ierr);
1087b77ba244SStefano Zampini     ierr = MatMult(shell->axpy,shell->axpy_right,shell->axpy_left);CHKERRQ(ierr);
1088b77ba244SStefano Zampini     ierr = VecAXPY(y,shell->axpy_vscale,shell->axpy_left);CHKERRQ(ierr);
10899f137db4SBarry Smith   }
1090ef66eb69SBarry Smith   PetscFunctionReturn(0);
1091ef66eb69SBarry Smith }
1092ef66eb69SBarry Smith 
10935edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
10945edf6516SJed Brown {
10955edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
10965edf6516SJed Brown   PetscErrorCode ierr;
10975edf6516SJed Brown 
10985edf6516SJed Brown   PetscFunctionBegin;
10995edf6516SJed Brown   if (y == z) {
11005edf6516SJed Brown     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
11015edf6516SJed Brown     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
1102b8430d49SBarry Smith     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
11035edf6516SJed Brown   } else {
11045edf6516SJed Brown     ierr = MatMult(A,x,z);CHKERRQ(ierr);
11055edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
11065edf6516SJed Brown   }
11075edf6516SJed Brown   PetscFunctionReturn(0);
11085edf6516SJed Brown }
11095edf6516SJed Brown 
111018be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
111118be62a5SSatish Balay {
111218be62a5SSatish Balay   Mat_Shell        *shell = (Mat_Shell*)A->data;
111318be62a5SSatish Balay   PetscErrorCode   ierr;
111425578ef6SJed Brown   Vec              xx;
1115e3f487b0SBarry Smith   PetscObjectState instate,outstate;
111618be62a5SSatish Balay 
111718be62a5SSatish Balay   PetscFunctionBegin;
1118b77ba244SStefano Zampini   if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL");
111945960306SStefano Zampini   ierr = MatShellPreZeroLeft(A,x,&xx);CHKERRQ(ierr);
112045960306SStefano Zampini   ierr = MatShellPreScaleLeft(A,xx,&xx);CHKERRQ(ierr);
1121e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
112228f357e3SAlex Fikl   ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr);
1123e3f487b0SBarry Smith   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
1124e3f487b0SBarry Smith   if (instate == outstate) {
1125e3f487b0SBarry Smith     /* increase the state of the output vector since the user did not update its state themself as should have been done */
1126e3f487b0SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
1127e3f487b0SBarry Smith   }
11288fe8eb27SJed Brown   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
11298fe8eb27SJed Brown   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
113045960306SStefano Zampini   ierr = MatShellPostZeroRight(A,y);CHKERRQ(ierr);
1131350ec83bSStefano Zampini 
1132350ec83bSStefano Zampini   if (shell->axpy) {
1133ef5c7bd2SStefano Zampini     Mat              X;
1134ef5c7bd2SStefano Zampini     PetscObjectState axpy_state;
1135ef5c7bd2SStefano Zampini 
1136ef5c7bd2SStefano Zampini     ierr = MatShellGetContext(shell->axpy,(void *)&X);CHKERRQ(ierr);
1137ef5c7bd2SStefano Zampini     ierr = PetscObjectStateGet((PetscObject)X,&axpy_state);CHKERRQ(ierr);
113831c70108SStefano Zampini     if (shell->axpy_state != axpy_state) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)");
1139b77ba244SStefano Zampini     ierr = MatCreateVecs(shell->axpy,shell->axpy_right ? NULL : &shell->axpy_right,shell->axpy_left ? NULL : &shell->axpy_left);CHKERRQ(ierr);
1140b77ba244SStefano Zampini     ierr = VecCopy(x,shell->axpy_left);CHKERRQ(ierr);
1141b77ba244SStefano Zampini     ierr = MatMultTranspose(shell->axpy,shell->axpy_left,shell->axpy_right);CHKERRQ(ierr);
1142b77ba244SStefano Zampini     ierr = VecAXPY(y,shell->axpy_vscale,shell->axpy_right);CHKERRQ(ierr);
1143350ec83bSStefano Zampini   }
114418be62a5SSatish Balay   PetscFunctionReturn(0);
114518be62a5SSatish Balay }
114618be62a5SSatish Balay 
11475edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
11485edf6516SJed Brown {
11495edf6516SJed Brown   Mat_Shell      *shell = (Mat_Shell*)A->data;
11505edf6516SJed Brown   PetscErrorCode ierr;
11515edf6516SJed Brown 
11525edf6516SJed Brown   PetscFunctionBegin;
11535edf6516SJed Brown   if (y == z) {
11545edf6516SJed Brown     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
11555edf6516SJed Brown     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
1156dac989eaS“Marek     ierr = VecAXPY(z,1.0,shell->left_add_work);CHKERRQ(ierr);
11575edf6516SJed Brown   } else {
11585edf6516SJed Brown     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
11595edf6516SJed Brown     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
11605edf6516SJed Brown   }
11615edf6516SJed Brown   PetscFunctionReturn(0);
11625edf6516SJed Brown }
11635edf6516SJed Brown 
11640c0fd78eSBarry Smith /*
11650c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
11660c0fd78eSBarry Smith */
116718be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
116818be62a5SSatish Balay {
116918be62a5SSatish Balay   Mat_Shell      *shell = (Mat_Shell*)A->data;
117018be62a5SSatish Balay   PetscErrorCode ierr;
117118be62a5SSatish Balay 
117218be62a5SSatish Balay   PetscFunctionBegin;
11730c0fd78eSBarry Smith   if (shell->ops->getdiagonal) {
117428f357e3SAlex Fikl     ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr);
1175305211d5SBarry 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,...)");
117618be62a5SSatish Balay   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
11778fe8eb27SJed Brown   if (shell->dshift) {
11780c0fd78eSBarry Smith     ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr);
117918be62a5SSatish Balay   }
11800c0fd78eSBarry Smith   ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
11818fe8eb27SJed Brown   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
11828fe8eb27SJed Brown   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
118345960306SStefano Zampini   if (shell->zrows) {
118445960306SStefano Zampini     ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
118545960306SStefano Zampini     ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
118645960306SStefano Zampini   }
118781c02519SBarry Smith   if (shell->axpy) {
1188ef5c7bd2SStefano Zampini     Mat              X;
1189ef5c7bd2SStefano Zampini     PetscObjectState axpy_state;
1190ef5c7bd2SStefano Zampini 
1191ef5c7bd2SStefano Zampini     ierr = MatShellGetContext(shell->axpy,(void *)&X);CHKERRQ(ierr);
1192ef5c7bd2SStefano Zampini     ierr = PetscObjectStateGet((PetscObject)X,&axpy_state);CHKERRQ(ierr);
119331c70108SStefano Zampini     if (shell->axpy_state != axpy_state) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)");
1194b77ba244SStefano Zampini     ierr = MatCreateVecs(shell->axpy,NULL,shell->axpy_left ? NULL : &shell->axpy_left);CHKERRQ(ierr);
1195b77ba244SStefano Zampini     ierr = MatGetDiagonal(shell->axpy,shell->axpy_left);CHKERRQ(ierr);
1196b77ba244SStefano Zampini     ierr = VecAXPY(v,shell->axpy_vscale,shell->axpy_left);CHKERRQ(ierr);
119781c02519SBarry Smith   }
119818be62a5SSatish Balay   PetscFunctionReturn(0);
119918be62a5SSatish Balay }
120018be62a5SSatish Balay 
1201f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
1202ef66eb69SBarry Smith {
1203ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
12048fe8eb27SJed Brown   PetscErrorCode ierr;
1205789d8953SBarry Smith   PetscBool      flg;
1206b24ad042SBarry Smith 
1207ef66eb69SBarry Smith   PetscFunctionBegin;
1208789d8953SBarry Smith   ierr = MatHasCongruentLayouts(Y,&flg);CHKERRQ(ierr);
1209789d8953SBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot shift shell matrix if it is not congruent");
12100c0fd78eSBarry Smith   if (shell->left || shell->right) {
12118fe8eb27SJed Brown     if (!shell->dshift) {
12120c0fd78eSBarry Smith       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
12130c0fd78eSBarry Smith       ierr = VecSet(shell->dshift,a);CHKERRQ(ierr);
12140c0fd78eSBarry Smith     } else {
12150c0fd78eSBarry Smith       if (shell->left)  {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
12160c0fd78eSBarry Smith       if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
12179f137db4SBarry Smith       ierr = VecShift(shell->dshift,a);CHKERRQ(ierr);
12180c0fd78eSBarry Smith     }
12198fe8eb27SJed Brown     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
12208fe8eb27SJed Brown     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
12218fe8eb27SJed Brown   } else shell->vshift += a;
122245960306SStefano Zampini   if (shell->zrows) {
122345960306SStefano Zampini     ierr = VecShift(shell->zvals,a);CHKERRQ(ierr);
122445960306SStefano Zampini   }
1225ef66eb69SBarry Smith   PetscFunctionReturn(0);
1226ef66eb69SBarry Smith }
1227ef66eb69SBarry Smith 
1228b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s)
12296464bf51SAlex Fikl {
12306464bf51SAlex Fikl   Mat_Shell      *shell = (Mat_Shell*)A->data;
12316464bf51SAlex Fikl   PetscErrorCode ierr;
12326464bf51SAlex Fikl 
12336464bf51SAlex Fikl   PetscFunctionBegin;
12340c0fd78eSBarry Smith   if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);}
12350c0fd78eSBarry Smith   if (shell->left || shell->right) {
123692fabd1fSBarry Smith     if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);}
12379f137db4SBarry Smith     if (shell->left && shell->right)  {
12389f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
12399f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
12409f137db4SBarry Smith     } else if (shell->left) {
12419f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
12429f137db4SBarry Smith     } else {
12439f137db4SBarry Smith       ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr);
12449f137db4SBarry Smith     }
1245b253ae0bS“Marek     ierr = VecAXPY(shell->dshift,s,shell->right_work);CHKERRQ(ierr);
12460c0fd78eSBarry Smith   } else {
1247b253ae0bS“Marek     ierr = VecAXPY(shell->dshift,s,D);CHKERRQ(ierr);
1248b253ae0bS“Marek   }
1249b253ae0bS“Marek   PetscFunctionReturn(0);
1250b253ae0bS“Marek }
1251b253ae0bS“Marek 
1252b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
1253b253ae0bS“Marek {
125445960306SStefano Zampini   Mat_Shell      *shell = (Mat_Shell*)A->data;
1255b253ae0bS“Marek   Vec            d;
1256b253ae0bS“Marek   PetscErrorCode ierr;
1257789d8953SBarry Smith   PetscBool      flg;
1258b253ae0bS“Marek 
1259b253ae0bS“Marek   PetscFunctionBegin;
1260789d8953SBarry Smith   ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr);
1261789d8953SBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot diagonal set or shift shell matrix if it is not congruent");
1262b253ae0bS“Marek   if (ins == INSERT_VALUES) {
1263b253ae0bS“Marek     if (!A->ops->getdiagonal) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first");
1264b253ae0bS“Marek     ierr = VecDuplicate(D,&d);CHKERRQ(ierr);
1265b253ae0bS“Marek     ierr = MatGetDiagonal(A,d);CHKERRQ(ierr);
1266b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,d,-1.);CHKERRQ(ierr);
1267b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr);
1268b253ae0bS“Marek     ierr = VecDestroy(&d);CHKERRQ(ierr);
126945960306SStefano Zampini     if (shell->zrows) {
127045960306SStefano Zampini       ierr = VecCopy(D,shell->zvals);CHKERRQ(ierr);
127145960306SStefano Zampini     }
1272b253ae0bS“Marek   } else {
1273b253ae0bS“Marek     ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr);
127445960306SStefano Zampini     if (shell->zrows) {
127545960306SStefano Zampini       ierr = VecAXPY(shell->zvals,1.0,D);CHKERRQ(ierr);
127645960306SStefano Zampini     }
12776464bf51SAlex Fikl   }
12786464bf51SAlex Fikl   PetscFunctionReturn(0);
12796464bf51SAlex Fikl }
12806464bf51SAlex Fikl 
1281f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
1282ef66eb69SBarry Smith {
1283ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
12848fe8eb27SJed Brown   PetscErrorCode ierr;
1285b24ad042SBarry Smith 
1286ef66eb69SBarry Smith   PetscFunctionBegin;
1287f4df32b1SMatthew Knepley   shell->vscale *= a;
12880c0fd78eSBarry Smith   shell->vshift *= a;
12892205254eSKarl Rupp   if (shell->dshift) {
12902205254eSKarl Rupp     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
12910c0fd78eSBarry Smith   }
129281c02519SBarry Smith   shell->axpy_vscale *= a;
129345960306SStefano Zampini   if (shell->zrows) {
129445960306SStefano Zampini     ierr = VecScale(shell->zvals,a);CHKERRQ(ierr);
129545960306SStefano Zampini   }
12968fe8eb27SJed Brown   PetscFunctionReturn(0);
129718be62a5SSatish Balay }
12988fe8eb27SJed Brown 
12998fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
13008fe8eb27SJed Brown {
13018fe8eb27SJed Brown   Mat_Shell      *shell = (Mat_Shell*)Y->data;
13028fe8eb27SJed Brown   PetscErrorCode ierr;
13038fe8eb27SJed Brown 
13048fe8eb27SJed Brown   PetscFunctionBegin;
13058fe8eb27SJed Brown   if (left) {
13060c0fd78eSBarry Smith     if (!shell->left) {
13070c0fd78eSBarry Smith       ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr);
13088fe8eb27SJed Brown       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
13090c0fd78eSBarry Smith     } else {
13100c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
131118be62a5SSatish Balay     }
131245960306SStefano Zampini     if (shell->zrows) {
131345960306SStefano Zampini       ierr = VecPointwiseMult(shell->zvals,shell->zvals,left);CHKERRQ(ierr);
131445960306SStefano Zampini     }
1315ef66eb69SBarry Smith   }
13168fe8eb27SJed Brown   if (right) {
13170c0fd78eSBarry Smith     if (!shell->right) {
13180c0fd78eSBarry Smith       ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr);
13198fe8eb27SJed Brown       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
13200c0fd78eSBarry Smith     } else {
13210c0fd78eSBarry Smith       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
13228fe8eb27SJed Brown     }
132345960306SStefano Zampini     if (shell->zrows) {
132445960306SStefano Zampini       if (!shell->left_work) {
132545960306SStefano Zampini         ierr = MatCreateVecs(Y,NULL,&shell->left_work);CHKERRQ(ierr);
132645960306SStefano Zampini       }
132745960306SStefano Zampini       ierr = VecSet(shell->zvals_w,1.0);CHKERRQ(ierr);
132845960306SStefano Zampini       ierr = VecScatterBegin(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
132945960306SStefano Zampini       ierr = VecScatterEnd(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
133045960306SStefano Zampini       ierr = VecPointwiseMult(shell->zvals,shell->zvals,shell->zvals_w);CHKERRQ(ierr);
133145960306SStefano Zampini     }
13328fe8eb27SJed Brown   }
1333ef5c7bd2SStefano Zampini   if (shell->axpy) {
1334ef5c7bd2SStefano Zampini     ierr = MatDiagonalScale(shell->axpy,left,right);CHKERRQ(ierr);
1335ef5c7bd2SStefano Zampini   }
1336ef66eb69SBarry Smith   PetscFunctionReturn(0);
1337ef66eb69SBarry Smith }
1338ef66eb69SBarry Smith 
1339dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
1340ef66eb69SBarry Smith {
1341ef66eb69SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
13420c0fd78eSBarry Smith   PetscErrorCode ierr;
1343ef66eb69SBarry Smith 
1344ef66eb69SBarry Smith   PetscFunctionBegin;
13458fe8eb27SJed Brown   if (t == MAT_FINAL_ASSEMBLY) {
1346ef66eb69SBarry Smith     shell->vshift = 0.0;
1347ef66eb69SBarry Smith     shell->vscale = 1.0;
1348ef5c7bd2SStefano Zampini     shell->axpy_vscale = 0.0;
1349ef5c7bd2SStefano Zampini     shell->axpy_state  = 0;
13500c0fd78eSBarry Smith     ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
13510c0fd78eSBarry Smith     ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
13520c0fd78eSBarry Smith     ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
135381c02519SBarry Smith     ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
1354b77ba244SStefano Zampini     ierr = VecDestroy(&shell->axpy_left);CHKERRQ(ierr);
1355b77ba244SStefano Zampini     ierr = VecDestroy(&shell->axpy_right);CHKERRQ(ierr);
135645960306SStefano Zampini     ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr);
135745960306SStefano Zampini     ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr);
135845960306SStefano Zampini     ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr);
135945960306SStefano Zampini     ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr);
1360ef66eb69SBarry Smith   }
1361ef66eb69SBarry Smith   PetscFunctionReturn(0);
1362ef66eb69SBarry Smith }
1363ef66eb69SBarry Smith 
13643b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d)
13653b49f96aSBarry Smith {
13663b49f96aSBarry Smith   PetscFunctionBegin;
13673b49f96aSBarry Smith   *missing = PETSC_FALSE;
13683b49f96aSBarry Smith   PetscFunctionReturn(0);
13693b49f96aSBarry Smith }
13703b49f96aSBarry Smith 
13719f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
13729f137db4SBarry Smith {
13739f137db4SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)Y->data;
13749f137db4SBarry Smith   PetscErrorCode ierr;
13759f137db4SBarry Smith 
13769f137db4SBarry Smith   PetscFunctionBegin;
1377ef5c7bd2SStefano Zampini   if (X == Y) {
1378ef5c7bd2SStefano Zampini     ierr = MatScale(Y,1.0 + a);CHKERRQ(ierr);
1379ef5c7bd2SStefano Zampini     PetscFunctionReturn(0);
1380ef5c7bd2SStefano Zampini   }
1381ef5c7bd2SStefano Zampini   if (!shell->axpy) {
1382ef5c7bd2SStefano Zampini     ierr = MatConvertFrom_Shell(X,MATSHELL,MAT_INITIAL_MATRIX,&shell->axpy);CHKERRQ(ierr);
13839f137db4SBarry Smith     shell->axpy_vscale = a;
1384ef5c7bd2SStefano Zampini     ierr = PetscObjectStateGet((PetscObject)X,&shell->axpy_state);CHKERRQ(ierr);
1385ef5c7bd2SStefano Zampini   } else {
1386ef5c7bd2SStefano Zampini     ierr = MatAXPY(shell->axpy,a/shell->axpy_vscale,X,str);CHKERRQ(ierr);
1387ef5c7bd2SStefano Zampini   }
13889f137db4SBarry Smith   PetscFunctionReturn(0);
13899f137db4SBarry Smith }
13909f137db4SBarry Smith 
1391f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL,
1392f4259b30SLisandro Dalcin                                        NULL,
1393f4259b30SLisandro Dalcin                                        NULL,
1394f4259b30SLisandro Dalcin                                        NULL,
13950c0fd78eSBarry Smith                                 /* 4*/ MatMultAdd_Shell,
1396f4259b30SLisandro Dalcin                                        NULL,
13970c0fd78eSBarry Smith                                        MatMultTransposeAdd_Shell,
1398f4259b30SLisandro Dalcin                                        NULL,
1399f4259b30SLisandro Dalcin                                        NULL,
1400f4259b30SLisandro Dalcin                                        NULL,
1401f4259b30SLisandro Dalcin                                 /*10*/ NULL,
1402f4259b30SLisandro Dalcin                                        NULL,
1403f4259b30SLisandro Dalcin                                        NULL,
1404f4259b30SLisandro Dalcin                                        NULL,
1405f4259b30SLisandro Dalcin                                        NULL,
1406f4259b30SLisandro Dalcin                                 /*15*/ NULL,
1407f4259b30SLisandro Dalcin                                        NULL,
1408f4259b30SLisandro Dalcin                                        NULL,
14098fe8eb27SJed Brown                                        MatDiagonalScale_Shell,
1410f4259b30SLisandro Dalcin                                        NULL,
1411f4259b30SLisandro Dalcin                                 /*20*/ NULL,
1412ef66eb69SBarry Smith                                        MatAssemblyEnd_Shell,
1413f4259b30SLisandro Dalcin                                        NULL,
1414f4259b30SLisandro Dalcin                                        NULL,
141545960306SStefano Zampini                                 /*24*/ MatZeroRows_Shell,
1416f4259b30SLisandro Dalcin                                        NULL,
1417f4259b30SLisandro Dalcin                                        NULL,
1418f4259b30SLisandro Dalcin                                        NULL,
1419f4259b30SLisandro Dalcin                                        NULL,
1420f4259b30SLisandro Dalcin                                 /*29*/ NULL,
1421f4259b30SLisandro Dalcin                                        NULL,
1422f4259b30SLisandro Dalcin                                        NULL,
1423f4259b30SLisandro Dalcin                                        NULL,
1424f4259b30SLisandro Dalcin                                        NULL,
1425cb8c8a77SBarry Smith                                 /*34*/ MatDuplicate_Shell,
1426f4259b30SLisandro Dalcin                                        NULL,
1427f4259b30SLisandro Dalcin                                        NULL,
1428f4259b30SLisandro Dalcin                                        NULL,
1429f4259b30SLisandro Dalcin                                        NULL,
14309f137db4SBarry Smith                                 /*39*/ MatAXPY_Shell,
1431f4259b30SLisandro Dalcin                                        NULL,
1432f4259b30SLisandro Dalcin                                        NULL,
1433f4259b30SLisandro Dalcin                                        NULL,
1434cb8c8a77SBarry Smith                                        MatCopy_Shell,
1435f4259b30SLisandro Dalcin                                 /*44*/ NULL,
1436ef66eb69SBarry Smith                                        MatScale_Shell,
1437ef66eb69SBarry Smith                                        MatShift_Shell,
14386464bf51SAlex Fikl                                        MatDiagonalSet_Shell,
143945960306SStefano Zampini                                        MatZeroRowsColumns_Shell,
1440f4259b30SLisandro Dalcin                                 /*49*/ NULL,
1441f4259b30SLisandro Dalcin                                        NULL,
1442f4259b30SLisandro Dalcin                                        NULL,
1443f4259b30SLisandro Dalcin                                        NULL,
1444f4259b30SLisandro Dalcin                                        NULL,
1445f4259b30SLisandro Dalcin                                 /*54*/ NULL,
1446f4259b30SLisandro Dalcin                                        NULL,
1447f4259b30SLisandro Dalcin                                        NULL,
1448f4259b30SLisandro Dalcin                                        NULL,
1449f4259b30SLisandro Dalcin                                        NULL,
1450f4259b30SLisandro Dalcin                                 /*59*/ NULL,
1451b9b97703SBarry Smith                                        MatDestroy_Shell,
1452f4259b30SLisandro Dalcin                                        NULL,
1453251fad3fSStefano Zampini                                        MatConvertFrom_Shell,
1454f4259b30SLisandro Dalcin                                        NULL,
1455f4259b30SLisandro Dalcin                                 /*64*/ NULL,
1456f4259b30SLisandro Dalcin                                        NULL,
1457f4259b30SLisandro Dalcin                                        NULL,
1458f4259b30SLisandro Dalcin                                        NULL,
1459f4259b30SLisandro Dalcin                                        NULL,
1460f4259b30SLisandro Dalcin                                 /*69*/ NULL,
1461f4259b30SLisandro Dalcin                                        NULL,
1462c87e5d42SMatthew Knepley                                        MatConvert_Shell,
1463f4259b30SLisandro Dalcin                                        NULL,
1464f4259b30SLisandro Dalcin                                        NULL,
1465f4259b30SLisandro Dalcin                                 /*74*/ NULL,
1466f4259b30SLisandro Dalcin                                        NULL,
1467f4259b30SLisandro Dalcin                                        NULL,
1468f4259b30SLisandro Dalcin                                        NULL,
1469f4259b30SLisandro Dalcin                                        NULL,
1470f4259b30SLisandro Dalcin                                 /*79*/ NULL,
1471f4259b30SLisandro Dalcin                                        NULL,
1472f4259b30SLisandro Dalcin                                        NULL,
1473f4259b30SLisandro Dalcin                                        NULL,
1474f4259b30SLisandro Dalcin                                        NULL,
1475f4259b30SLisandro Dalcin                                 /*84*/ NULL,
1476f4259b30SLisandro Dalcin                                        NULL,
1477f4259b30SLisandro Dalcin                                        NULL,
1478f4259b30SLisandro Dalcin                                        NULL,
1479f4259b30SLisandro Dalcin                                        NULL,
1480f4259b30SLisandro Dalcin                                 /*89*/ NULL,
1481f4259b30SLisandro Dalcin                                        NULL,
1482f4259b30SLisandro Dalcin                                        NULL,
1483f4259b30SLisandro Dalcin                                        NULL,
1484f4259b30SLisandro Dalcin                                        NULL,
1485f4259b30SLisandro Dalcin                                 /*94*/ NULL,
1486f4259b30SLisandro Dalcin                                        NULL,
1487f4259b30SLisandro Dalcin                                        NULL,
1488f4259b30SLisandro Dalcin                                        NULL,
1489f4259b30SLisandro Dalcin                                        NULL,
1490f4259b30SLisandro Dalcin                                 /*99*/ NULL,
1491f4259b30SLisandro Dalcin                                        NULL,
1492f4259b30SLisandro Dalcin                                        NULL,
1493f4259b30SLisandro Dalcin                                        NULL,
1494f4259b30SLisandro Dalcin                                        NULL,
1495f4259b30SLisandro Dalcin                                /*104*/ NULL,
1496f4259b30SLisandro Dalcin                                        NULL,
1497f4259b30SLisandro Dalcin                                        NULL,
1498f4259b30SLisandro Dalcin                                        NULL,
1499f4259b30SLisandro Dalcin                                        NULL,
1500f4259b30SLisandro Dalcin                                /*109*/ NULL,
1501f4259b30SLisandro Dalcin                                        NULL,
1502f4259b30SLisandro Dalcin                                        NULL,
1503f4259b30SLisandro Dalcin                                        NULL,
15043b49f96aSBarry Smith                                        MatMissingDiagonal_Shell,
1505f4259b30SLisandro Dalcin                                /*114*/ NULL,
1506f4259b30SLisandro Dalcin                                        NULL,
1507f4259b30SLisandro Dalcin                                        NULL,
1508f4259b30SLisandro Dalcin                                        NULL,
1509f4259b30SLisandro Dalcin                                        NULL,
1510f4259b30SLisandro Dalcin                                /*119*/ NULL,
1511f4259b30SLisandro Dalcin                                        NULL,
1512f4259b30SLisandro Dalcin                                        NULL,
1513f4259b30SLisandro Dalcin                                        NULL,
1514f4259b30SLisandro Dalcin                                        NULL,
1515f4259b30SLisandro Dalcin                                /*124*/ NULL,
1516f4259b30SLisandro Dalcin                                        NULL,
1517f4259b30SLisandro Dalcin                                        NULL,
1518f4259b30SLisandro Dalcin                                        NULL,
1519f4259b30SLisandro Dalcin                                        NULL,
1520f4259b30SLisandro Dalcin                                /*129*/ NULL,
1521f4259b30SLisandro Dalcin                                        NULL,
1522f4259b30SLisandro Dalcin                                        NULL,
1523f4259b30SLisandro Dalcin                                        NULL,
1524f4259b30SLisandro Dalcin                                        NULL,
1525f4259b30SLisandro Dalcin                                /*134*/ NULL,
1526f4259b30SLisandro Dalcin                                        NULL,
1527f4259b30SLisandro Dalcin                                        NULL,
1528f4259b30SLisandro Dalcin                                        NULL,
1529f4259b30SLisandro Dalcin                                        NULL,
1530f4259b30SLisandro Dalcin                                /*139*/ NULL,
1531f4259b30SLisandro Dalcin                                        NULL,
1532f4259b30SLisandro Dalcin                                        NULL
15333964eb88SJed Brown };
1534273d9f13SBarry Smith 
1535789d8953SBarry Smith PetscErrorCode  MatShellSetContext_Shell(Mat mat,void *ctx)
1536789d8953SBarry Smith {
1537789d8953SBarry Smith   Mat_Shell *shell = (Mat_Shell*)mat->data;
1538789d8953SBarry Smith 
1539789d8953SBarry Smith   PetscFunctionBegin;
1540789d8953SBarry Smith   shell->ctx = ctx;
1541789d8953SBarry Smith   PetscFunctionReturn(0);
1542789d8953SBarry Smith }
1543789d8953SBarry Smith 
1544db77db73SJeremy L Thompson static PetscErrorCode MatShellSetVecType_Shell(Mat mat,VecType vtype)
1545db77db73SJeremy L Thompson {
1546db77db73SJeremy L Thompson   PetscErrorCode ierr;
1547db77db73SJeremy L Thompson 
1548db77db73SJeremy L Thompson   PetscFunctionBegin;
1549db77db73SJeremy L Thompson   ierr = PetscFree(mat->defaultvectype);CHKERRQ(ierr);
1550db77db73SJeremy L Thompson   ierr = PetscStrallocpy(vtype,(char**)&mat->defaultvectype);CHKERRQ(ierr);
1551db77db73SJeremy L Thompson   PetscFunctionReturn(0);
1552db77db73SJeremy L Thompson }
1553db77db73SJeremy L Thompson 
1554789d8953SBarry Smith PetscErrorCode MatShellSetManageScalingShifts_Shell(Mat A)
1555789d8953SBarry Smith {
1556789d8953SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)A->data;
1557789d8953SBarry Smith 
1558789d8953SBarry Smith   PetscFunctionBegin;
1559789d8953SBarry Smith   shell->managescalingshifts = PETSC_FALSE;
1560789d8953SBarry Smith   A->ops->diagonalset   = NULL;
1561789d8953SBarry Smith   A->ops->diagonalscale = NULL;
1562789d8953SBarry Smith   A->ops->scale         = NULL;
1563789d8953SBarry Smith   A->ops->shift         = NULL;
1564789d8953SBarry Smith   A->ops->axpy          = NULL;
1565789d8953SBarry Smith   PetscFunctionReturn(0);
1566789d8953SBarry Smith }
1567789d8953SBarry Smith 
1568789d8953SBarry Smith PetscErrorCode MatShellSetOperation_Shell(Mat mat,MatOperation op,void (*f)(void))
1569789d8953SBarry Smith {
1570feb237baSPierre Jolivet   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1571789d8953SBarry Smith 
1572789d8953SBarry Smith   PetscFunctionBegin;
1573789d8953SBarry Smith   switch (op) {
1574789d8953SBarry Smith   case MATOP_DESTROY:
1575789d8953SBarry Smith     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
1576789d8953SBarry Smith     break;
1577789d8953SBarry Smith   case MATOP_VIEW:
1578789d8953SBarry Smith     if (!mat->ops->viewnative) {
1579789d8953SBarry Smith       mat->ops->viewnative = mat->ops->view;
1580789d8953SBarry Smith     }
1581789d8953SBarry Smith     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
1582789d8953SBarry Smith     break;
1583789d8953SBarry Smith   case MATOP_COPY:
1584789d8953SBarry Smith     shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f;
1585789d8953SBarry Smith     break;
1586789d8953SBarry Smith   case MATOP_DIAGONAL_SET:
1587789d8953SBarry Smith   case MATOP_DIAGONAL_SCALE:
1588789d8953SBarry Smith   case MATOP_SHIFT:
1589789d8953SBarry Smith   case MATOP_SCALE:
1590789d8953SBarry Smith   case MATOP_AXPY:
1591789d8953SBarry Smith   case MATOP_ZERO_ROWS:
1592789d8953SBarry Smith   case MATOP_ZERO_ROWS_COLUMNS:
1593789d8953SBarry Smith     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
1594789d8953SBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
1595789d8953SBarry Smith     break;
1596789d8953SBarry Smith   case MATOP_GET_DIAGONAL:
1597789d8953SBarry Smith     if (shell->managescalingshifts) {
1598789d8953SBarry Smith       shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1599789d8953SBarry Smith       mat->ops->getdiagonal   = MatGetDiagonal_Shell;
1600789d8953SBarry Smith     } else {
1601789d8953SBarry Smith       shell->ops->getdiagonal = NULL;
1602789d8953SBarry Smith       mat->ops->getdiagonal   = (PetscErrorCode (*)(Mat,Vec))f;
1603789d8953SBarry Smith     }
1604789d8953SBarry Smith     break;
1605789d8953SBarry Smith   case MATOP_MULT:
1606789d8953SBarry Smith     if (shell->managescalingshifts) {
1607789d8953SBarry Smith       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1608789d8953SBarry Smith       mat->ops->mult   = MatMult_Shell;
1609789d8953SBarry Smith     } else {
1610789d8953SBarry Smith       shell->ops->mult = NULL;
1611789d8953SBarry Smith       mat->ops->mult   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1612789d8953SBarry Smith     }
1613789d8953SBarry Smith     break;
1614789d8953SBarry Smith   case MATOP_MULT_TRANSPOSE:
1615789d8953SBarry Smith     if (shell->managescalingshifts) {
1616789d8953SBarry Smith       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1617789d8953SBarry Smith       mat->ops->multtranspose   = MatMultTranspose_Shell;
1618789d8953SBarry Smith     } else {
1619789d8953SBarry Smith       shell->ops->multtranspose = NULL;
1620789d8953SBarry Smith       mat->ops->multtranspose   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1621789d8953SBarry Smith     }
1622789d8953SBarry Smith     break;
1623789d8953SBarry Smith   default:
1624789d8953SBarry Smith     (((void(**)(void))mat->ops)[op]) = f;
1625789d8953SBarry Smith     break;
1626789d8953SBarry Smith   }
1627789d8953SBarry Smith   PetscFunctionReturn(0);
1628789d8953SBarry Smith }
1629789d8953SBarry Smith 
1630789d8953SBarry Smith PetscErrorCode MatShellGetOperation_Shell(Mat mat,MatOperation op,void(**f)(void))
1631789d8953SBarry Smith {
1632789d8953SBarry Smith   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1633789d8953SBarry Smith 
1634789d8953SBarry Smith   PetscFunctionBegin;
1635789d8953SBarry Smith   switch (op) {
1636789d8953SBarry Smith   case MATOP_DESTROY:
1637789d8953SBarry Smith     *f = (void (*)(void))shell->ops->destroy;
1638789d8953SBarry Smith     break;
1639789d8953SBarry Smith   case MATOP_VIEW:
1640789d8953SBarry Smith     *f = (void (*)(void))mat->ops->view;
1641789d8953SBarry Smith     break;
1642789d8953SBarry Smith   case MATOP_COPY:
1643789d8953SBarry Smith     *f = (void (*)(void))shell->ops->copy;
1644789d8953SBarry Smith     break;
1645789d8953SBarry Smith   case MATOP_DIAGONAL_SET:
1646789d8953SBarry Smith   case MATOP_DIAGONAL_SCALE:
1647789d8953SBarry Smith   case MATOP_SHIFT:
1648789d8953SBarry Smith   case MATOP_SCALE:
1649789d8953SBarry Smith   case MATOP_AXPY:
1650789d8953SBarry Smith   case MATOP_ZERO_ROWS:
1651789d8953SBarry Smith   case MATOP_ZERO_ROWS_COLUMNS:
1652789d8953SBarry Smith     *f = (((void (**)(void))mat->ops)[op]);
1653789d8953SBarry Smith     break;
1654789d8953SBarry Smith   case MATOP_GET_DIAGONAL:
1655789d8953SBarry Smith     if (shell->ops->getdiagonal)
1656789d8953SBarry Smith       *f = (void (*)(void))shell->ops->getdiagonal;
1657789d8953SBarry Smith     else
1658789d8953SBarry Smith       *f = (((void (**)(void))mat->ops)[op]);
1659789d8953SBarry Smith     break;
1660789d8953SBarry Smith   case MATOP_MULT:
1661789d8953SBarry Smith     if (shell->ops->mult)
1662789d8953SBarry Smith       *f = (void (*)(void))shell->ops->mult;
1663789d8953SBarry Smith     else
1664789d8953SBarry Smith       *f = (((void (**)(void))mat->ops)[op]);
1665789d8953SBarry Smith     break;
1666789d8953SBarry Smith   case MATOP_MULT_TRANSPOSE:
1667789d8953SBarry Smith     if (shell->ops->multtranspose)
1668789d8953SBarry Smith       *f = (void (*)(void))shell->ops->multtranspose;
1669789d8953SBarry Smith     else
1670789d8953SBarry Smith       *f = (((void (**)(void))mat->ops)[op]);
1671789d8953SBarry Smith     break;
1672789d8953SBarry Smith   default:
1673789d8953SBarry Smith     *f = (((void (**)(void))mat->ops)[op]);
1674789d8953SBarry Smith   }
1675789d8953SBarry Smith   PetscFunctionReturn(0);
1676789d8953SBarry Smith }
1677789d8953SBarry Smith 
16780bad9183SKris Buschelman /*MC
1679fafad747SKris Buschelman    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
16800bad9183SKris Buschelman 
16810bad9183SKris Buschelman   Level: advanced
16820bad9183SKris Buschelman 
16830c0fd78eSBarry Smith .seealso: MatCreateShell()
16840bad9183SKris Buschelman M*/
16850bad9183SKris Buschelman 
16868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
1687273d9f13SBarry Smith {
1688273d9f13SBarry Smith   Mat_Shell      *b;
1689dfbe8321SBarry Smith   PetscErrorCode ierr;
1690273d9f13SBarry Smith 
1691273d9f13SBarry Smith   PetscFunctionBegin;
1692273d9f13SBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1693273d9f13SBarry Smith 
1694b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1695273d9f13SBarry Smith   A->data = (void*)b;
1696273d9f13SBarry Smith 
1697f4259b30SLisandro Dalcin   b->ctx                 = NULL;
1698ef66eb69SBarry Smith   b->vshift              = 0.0;
1699ef66eb69SBarry Smith   b->vscale              = 1.0;
17000c0fd78eSBarry Smith   b->managescalingshifts = PETSC_TRUE;
1701273d9f13SBarry Smith   A->assembled           = PETSC_TRUE;
1702210f0121SBarry Smith   A->preallocated        = PETSC_FALSE;
17032205254eSKarl Rupp 
1704789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetContext_C",MatShellGetContext_Shell);CHKERRQ(ierr);
1705789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetContext_C",MatShellSetContext_Shell);CHKERRQ(ierr);
1706db77db73SJeremy L Thompson   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetVecType_C",MatShellSetVecType_Shell);CHKERRQ(ierr);
1707789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetManageScalingShifts_C",MatShellSetManageScalingShifts_Shell);CHKERRQ(ierr);
1708789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetOperation_C",MatShellSetOperation_Shell);CHKERRQ(ierr);
1709789d8953SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetOperation_C",MatShellGetOperation_Shell);CHKERRQ(ierr);
1710b77ba244SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetMatProductOperation_C",MatShellSetMatProductOperation_Shell);CHKERRQ(ierr);
171117667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
1712273d9f13SBarry Smith   PetscFunctionReturn(0);
1713273d9f13SBarry Smith }
1714e51e0e81SBarry Smith 
17154b828684SBarry Smith /*@C
1716052efed2SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
1717ff756334SLois Curfman McInnes    private data storage format.
1718e51e0e81SBarry Smith 
1719d083f849SBarry Smith   Collective
1720c7fcc2eaSBarry Smith 
1721e51e0e81SBarry Smith    Input Parameters:
1722c7fcc2eaSBarry Smith +  comm - MPI communicator
1723c7fcc2eaSBarry Smith .  m - number of local rows (must be given)
1724c7fcc2eaSBarry Smith .  n - number of local columns (must be given)
1725c7fcc2eaSBarry Smith .  M - number of global rows (may be PETSC_DETERMINE)
1726c7fcc2eaSBarry Smith .  N - number of global columns (may be PETSC_DETERMINE)
1727c7fcc2eaSBarry Smith -  ctx - pointer to data needed by the shell matrix routines
1728e51e0e81SBarry Smith 
1729ff756334SLois Curfman McInnes    Output Parameter:
173044cd7ae7SLois Curfman McInnes .  A - the matrix
1731e51e0e81SBarry Smith 
1732ff2fd236SBarry Smith    Level: advanced
1733ff2fd236SBarry Smith 
1734f39d1f56SLois Curfman McInnes   Usage:
17355bd1e576SStefano Zampini $    extern PetscErrorCode mult(Mat,Vec,Vec);
1736f39d1f56SLois Curfman McInnes $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
1737c134de8dSSatish Balay $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
1738f39d1f56SLois Curfman McInnes $    [ Use matrix for operations that have been set ]
1739f39d1f56SLois Curfman McInnes $    MatDestroy(mat);
1740f39d1f56SLois Curfman McInnes 
1741ff756334SLois Curfman McInnes    Notes:
1742ff756334SLois Curfman McInnes    The shell matrix type is intended to provide a simple class to use
1743ff756334SLois Curfman McInnes    with KSP (such as, for use with matrix-free methods). You should not
1744ff756334SLois Curfman McInnes    use the shell type if you plan to define a complete matrix class.
1745e51e0e81SBarry Smith 
174695452b02SPatrick Sanan    Fortran Notes:
174795452b02SPatrick Sanan     To use this from Fortran with a ctx you must write an interface definition for this
1748daf670e6SBarry Smith     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
1749daf670e6SBarry Smith     in as the ctx argument.
1750f38a8d56SBarry Smith 
1751f39d1f56SLois Curfman McInnes    PETSc requires that matrices and vectors being used for certain
1752f39d1f56SLois Curfman McInnes    operations are partitioned accordingly.  For example, when
1753645985a0SLois Curfman McInnes    creating a shell matrix, A, that supports parallel matrix-vector
1754645985a0SLois Curfman McInnes    products using MatMult(A,x,y) the user should set the number
1755645985a0SLois Curfman McInnes    of local matrix rows to be the number of local elements of the
1756645985a0SLois Curfman McInnes    corresponding result vector, y. Note that this is information is
1757645985a0SLois Curfman McInnes    required for use of the matrix interface routines, even though
1758645985a0SLois Curfman McInnes    the shell matrix may not actually be physically partitioned.
1759645985a0SLois Curfman McInnes    For example,
1760f39d1f56SLois Curfman McInnes 
1761f39d1f56SLois Curfman McInnes $
1762f39d1f56SLois Curfman McInnes $     Vec x, y
17635bd1e576SStefano Zampini $     extern PetscErrorCode mult(Mat,Vec,Vec);
1764645985a0SLois Curfman McInnes $     Mat A
1765f39d1f56SLois Curfman McInnes $
1766c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
1767c94f878dSBarry Smith $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
1768f39d1f56SLois Curfman McInnes $     VecGetLocalSize(y,&m);
1769c7fcc2eaSBarry Smith $     VecGetLocalSize(x,&n);
1770c7fcc2eaSBarry Smith $     MatCreateShell(comm,m,n,M,N,ctx,&A);
1771c134de8dSSatish Balay $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
1772645985a0SLois Curfman McInnes $     MatMult(A,x,y);
177345960306SStefano Zampini $     MatDestroy(&A);
177445960306SStefano Zampini $     VecDestroy(&y);
177545960306SStefano Zampini $     VecDestroy(&x);
1776645985a0SLois Curfman McInnes $
1777e51e0e81SBarry Smith 
17789b53d762SBarry Smith 
177945960306SStefano Zampini    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), MatScale(), MatZeroRows() and MatZeroRowsColumns() internally, so these
17809b53d762SBarry Smith    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.
17819b53d762SBarry Smith 
17829b53d762SBarry Smith 
17830c0fd78eSBarry Smith     For rectangular matrices do all the scalings and shifts make sense?
17840c0fd78eSBarry Smith 
178595452b02SPatrick Sanan     Developers Notes:
178695452b02SPatrick Sanan     Regarding shifting and scaling. The general form is
17870c0fd78eSBarry Smith 
17880c0fd78eSBarry Smith           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
17890c0fd78eSBarry Smith 
17900c0fd78eSBarry Smith       The order you apply the operations is important. For example if you have a dshift then
17910c0fd78eSBarry Smith       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
17920c0fd78eSBarry Smith       you get s*vscale*A + diag(shift)
17930c0fd78eSBarry Smith 
17940c0fd78eSBarry Smith           A is the user provided function.
17950c0fd78eSBarry Smith 
1796ad2f5c3fSBarry Smith    KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for
1797ad2f5c3fSBarry Smith    for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger
1798ad2f5c3fSBarry Smith    an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat);
1799ad2f5c3fSBarry Smith    each time the MATSHELL matrix has changed.
1800ad2f5c3fSBarry Smith 
1801b77ba244SStefano Zampini    Matrix product operations (i.e. MatMat, MatTranposeMat etc) can be specified using MatShellSetMatProductOperation()
1802b77ba244SStefano Zampini 
1803ad2f5c3fSBarry Smith    Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided
1804ad2f5c3fSBarry Smith    with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale().
1805ad2f5c3fSBarry Smith 
1806b77ba244SStefano Zampini .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts(), MatShellSetMatProductOperation()
1807e51e0e81SBarry Smith @*/
18087087cfbeSBarry Smith PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
1809e51e0e81SBarry Smith {
1810dfbe8321SBarry Smith   PetscErrorCode ierr;
1811ed3cc1f0SBarry Smith 
18123a40ed3dSBarry Smith   PetscFunctionBegin;
1813f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1814f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1815273d9f13SBarry Smith   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
1816273d9f13SBarry Smith   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
18170fcd0fccSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
1818273d9f13SBarry Smith   PetscFunctionReturn(0);
1819c7fcc2eaSBarry Smith }
1820c7fcc2eaSBarry Smith 
1821c6866cfdSSatish Balay /*@
1822273d9f13SBarry Smith     MatShellSetContext - sets the context for a shell matrix
1823c7fcc2eaSBarry Smith 
18243f9fe445SBarry Smith    Logically Collective on Mat
1825c7fcc2eaSBarry Smith 
1826273d9f13SBarry Smith     Input Parameters:
1827273d9f13SBarry Smith +   mat - the shell matrix
1828273d9f13SBarry Smith -   ctx - the context
1829273d9f13SBarry Smith 
1830273d9f13SBarry Smith    Level: advanced
1831273d9f13SBarry Smith 
183295452b02SPatrick Sanan    Fortran Notes:
183395452b02SPatrick Sanan     To use this from Fortran you must write a Fortran interface definition for this
1834daf670e6SBarry Smith     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1835273d9f13SBarry Smith 
1836273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
18370bc0a719SBarry Smith @*/
18387087cfbeSBarry Smith PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
1839273d9f13SBarry Smith {
1840dfbe8321SBarry Smith   PetscErrorCode ierr;
1841273d9f13SBarry Smith 
1842273d9f13SBarry Smith   PetscFunctionBegin;
18430700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1844ef5c7bd2SStefano Zampini   ierr = PetscTryMethod(mat,"MatShellSetContext_C",(Mat,void*),(mat,ctx));CHKERRQ(ierr);
18453a40ed3dSBarry Smith   PetscFunctionReturn(0);
1846e51e0e81SBarry Smith }
1847e51e0e81SBarry Smith 
1848db77db73SJeremy L Thompson /*@C
1849db77db73SJeremy L Thompson  MatShellSetVecType - Sets the type of Vec returned by MatCreateVecs()
1850db77db73SJeremy L Thompson 
1851db77db73SJeremy L Thompson  Logically collective
1852db77db73SJeremy L Thompson 
1853db77db73SJeremy L Thompson     Input Parameters:
1854db77db73SJeremy L Thompson +   mat   - the shell matrix
1855db77db73SJeremy L Thompson -   vtype - type to use for creating vectors
1856db77db73SJeremy L Thompson 
1857db77db73SJeremy L Thompson  Notes:
1858db77db73SJeremy L Thompson 
1859db77db73SJeremy L Thompson  Level: advanced
1860db77db73SJeremy L Thompson 
1861db77db73SJeremy L Thompson .seealso: MatCreateVecs()
1862db77db73SJeremy L Thompson @*/
1863db77db73SJeremy L Thompson PetscErrorCode  MatShellSetVecType(Mat mat,VecType vtype)
1864db77db73SJeremy L Thompson {
1865db77db73SJeremy L Thompson   PetscErrorCode ierr;
1866db77db73SJeremy L Thompson 
1867db77db73SJeremy L Thompson   PetscFunctionBegin;
1868db77db73SJeremy L Thompson   ierr = PetscTryMethod(mat,"MatShellSetVecType_C",(Mat,VecType),(mat,vtype));CHKERRQ(ierr);
1869db77db73SJeremy L Thompson   PetscFunctionReturn(0);
1870db77db73SJeremy L Thompson }
1871db77db73SJeremy L Thompson 
18720c0fd78eSBarry Smith /*@
18730c0fd78eSBarry Smith     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
18740c0fd78eSBarry Smith           after MatCreateShell()
18750c0fd78eSBarry Smith 
18760c0fd78eSBarry Smith    Logically Collective on Mat
18770c0fd78eSBarry Smith 
18780c0fd78eSBarry Smith     Input Parameter:
18790c0fd78eSBarry Smith .   mat - the shell matrix
18800c0fd78eSBarry Smith 
18810c0fd78eSBarry Smith   Level: advanced
18820c0fd78eSBarry Smith 
18830c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
18840c0fd78eSBarry Smith @*/
18850c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A)
18860c0fd78eSBarry Smith {
18870c0fd78eSBarry Smith   PetscErrorCode ierr;
18880c0fd78eSBarry Smith 
18890c0fd78eSBarry Smith   PetscFunctionBegin;
18900c0fd78eSBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1891ef5c7bd2SStefano Zampini   ierr = PetscTryMethod(A,"MatShellSetManageScalingShifts_C",(Mat),(A));CHKERRQ(ierr);
18920c0fd78eSBarry Smith   PetscFunctionReturn(0);
18930c0fd78eSBarry Smith }
18940c0fd78eSBarry Smith 
1895c16cb8f2SBarry Smith /*@C
1896f3b1f45cSBarry Smith     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
1897f3b1f45cSBarry Smith 
1898f3b1f45cSBarry Smith    Logically Collective on Mat
1899f3b1f45cSBarry Smith 
1900f3b1f45cSBarry Smith     Input Parameters:
1901f3b1f45cSBarry Smith +   mat - the shell matrix
1902f3b1f45cSBarry Smith .   f - the function
1903f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
1904f3b1f45cSBarry Smith -   ctx - an optional context for the function
1905f3b1f45cSBarry Smith 
1906f3b1f45cSBarry Smith    Output Parameter:
1907f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
1908f3b1f45cSBarry Smith 
1909f3b1f45cSBarry Smith    Options Database:
1910f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
1911f3b1f45cSBarry Smith 
1912f3b1f45cSBarry Smith    Level: advanced
1913f3b1f45cSBarry Smith 
191495452b02SPatrick Sanan    Fortran Notes:
191595452b02SPatrick Sanan     Not supported from Fortran
1916f3b1f45cSBarry Smith 
1917f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
1918f3b1f45cSBarry Smith @*/
1919f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
1920f3b1f45cSBarry Smith {
1921f3b1f45cSBarry Smith   PetscErrorCode ierr;
1922f3b1f45cSBarry Smith   PetscInt       m,n;
1923f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
1924f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
192574e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
1926f3b1f45cSBarry Smith 
1927f3b1f45cSBarry Smith   PetscFunctionBegin;
1928f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1929f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr);
1930f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
1931f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr);
1932f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
1933f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
1934f3b1f45cSBarry Smith 
19350bacdadaSStefano Zampini   ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr);
19360bacdadaSStefano Zampini   ierr = MatComputeOperator(mat,MATAIJ,&Dmat);CHKERRQ(ierr);
1937f3b1f45cSBarry Smith 
1938f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
1939f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
1940f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
1941f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
1942f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
1943f3b1f45cSBarry Smith     flag = PETSC_FALSE;
1944f3b1f45cSBarry Smith     if (v) {
1945fc7aafd1SBarry 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);
1946f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
1947f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
1948f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
1949f3b1f45cSBarry Smith     }
1950f3b1f45cSBarry Smith   } else if (v) {
1951fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
1952f3b1f45cSBarry Smith   }
1953f3b1f45cSBarry Smith   if (flg) *flg = flag;
1954f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
1955f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
1956f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
1957f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
1958f3b1f45cSBarry Smith   PetscFunctionReturn(0);
1959f3b1f45cSBarry Smith }
1960f3b1f45cSBarry Smith 
1961f3b1f45cSBarry Smith /*@C
1962f3b1f45cSBarry Smith     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.
1963f3b1f45cSBarry Smith 
1964f3b1f45cSBarry Smith    Logically Collective on Mat
1965f3b1f45cSBarry Smith 
1966f3b1f45cSBarry Smith     Input Parameters:
1967f3b1f45cSBarry Smith +   mat - the shell matrix
1968f3b1f45cSBarry Smith .   f - the function
1969f3b1f45cSBarry Smith .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
1970f3b1f45cSBarry Smith -   ctx - an optional context for the function
1971f3b1f45cSBarry Smith 
1972f3b1f45cSBarry Smith    Output Parameter:
1973f3b1f45cSBarry Smith .   flg - PETSC_TRUE if the multiply is likely correct
1974f3b1f45cSBarry Smith 
1975f3b1f45cSBarry Smith    Options Database:
1976f3b1f45cSBarry Smith .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
1977f3b1f45cSBarry Smith 
1978f3b1f45cSBarry Smith    Level: advanced
1979f3b1f45cSBarry Smith 
198095452b02SPatrick Sanan    Fortran Notes:
198195452b02SPatrick Sanan     Not supported from Fortran
1982f3b1f45cSBarry Smith 
1983f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
1984f3b1f45cSBarry Smith @*/
1985f3b1f45cSBarry Smith PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
1986f3b1f45cSBarry Smith {
1987f3b1f45cSBarry Smith   PetscErrorCode ierr;
1988f3b1f45cSBarry Smith   Vec            x,y,z;
1989f3b1f45cSBarry Smith   PetscInt       m,n,M,N;
1990f3b1f45cSBarry Smith   Mat            mf,Dmf,Dmat,Ddiff;
1991f3b1f45cSBarry Smith   PetscReal      Diffnorm,Dmfnorm;
199274e5cdcaSBarry Smith   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
1993f3b1f45cSBarry Smith 
1994f3b1f45cSBarry Smith   PetscFunctionBegin;
1995f3b1f45cSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1996f3b1f45cSBarry Smith   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
1997f3b1f45cSBarry Smith   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
1998f3b1f45cSBarry Smith   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
1999f3b1f45cSBarry Smith   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
2000f3b1f45cSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
2001f3b1f45cSBarry Smith   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
2002f3b1f45cSBarry Smith   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
2003f3b1f45cSBarry Smith   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
20040bacdadaSStefano Zampini   ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr);
2005f3b1f45cSBarry Smith   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
20060bacdadaSStefano Zampini   ierr = MatComputeOperatorTranspose(mat,MATAIJ,&Dmat);CHKERRQ(ierr);
2007f3b1f45cSBarry Smith 
2008f3b1f45cSBarry Smith   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
2009f3b1f45cSBarry Smith   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2010f3b1f45cSBarry Smith   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
2011f3b1f45cSBarry Smith   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
2012f3b1f45cSBarry Smith   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
2013f3b1f45cSBarry Smith     flag = PETSC_FALSE;
2014f3b1f45cSBarry Smith     if (v) {
2015fc7aafd1SBarry 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);
2016f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
2017f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
2018f3b1f45cSBarry Smith       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
2019f3b1f45cSBarry Smith     }
2020f3b1f45cSBarry Smith   } else if (v) {
2021fc7aafd1SBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
2022f3b1f45cSBarry Smith   }
2023f3b1f45cSBarry Smith   if (flg) *flg = flag;
2024f3b1f45cSBarry Smith   ierr = MatDestroy(&mf);CHKERRQ(ierr);
2025f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
2026f3b1f45cSBarry Smith   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
2027f3b1f45cSBarry Smith   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
2028f3b1f45cSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
2029f3b1f45cSBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
2030f3b1f45cSBarry Smith   ierr = VecDestroy(&z);CHKERRQ(ierr);
2031f3b1f45cSBarry Smith   PetscFunctionReturn(0);
2032f3b1f45cSBarry Smith }
2033f3b1f45cSBarry Smith 
2034f3b1f45cSBarry Smith /*@C
20350c0fd78eSBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
2036e51e0e81SBarry Smith 
20373f9fe445SBarry Smith    Logically Collective on Mat
2038fee21e36SBarry Smith 
2039c7fcc2eaSBarry Smith     Input Parameters:
2040c7fcc2eaSBarry Smith +   mat - the shell matrix
2041c7fcc2eaSBarry Smith .   op - the name of the operation
2042789d8953SBarry Smith -   g - the function that provides the operation.
2043c7fcc2eaSBarry Smith 
204415091d37SBarry Smith    Level: advanced
204515091d37SBarry Smith 
2046fae171e0SBarry Smith     Usage:
2047e93bc3c1Svictor $      extern PetscErrorCode usermult(Mat,Vec,Vec);
2048b77ba244SStefano Zampini $      MatCreateShell(comm,m,n,M,N,ctx,&A);
2049b77ba244SStefano Zampini $      MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
20500b627109SLois Curfman McInnes 
2051a62d957aSLois Curfman McInnes     Notes:
2052e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
20531c1c02c0SLois Curfman McInnes     operations, which all have the form MATOP_<OPERATION>, where
2054a62d957aSLois Curfman McInnes     <OPERATION> is the name (in all capital letters) of the
20551c1c02c0SLois Curfman McInnes     user interface routine (e.g., MatMult() -> MATOP_MULT).
2056a62d957aSLois Curfman McInnes 
2057a4d85fd7SBarry Smith     All user-provided functions (except for MATOP_DESTROY) should have the same calling
2058deebb3c3SLois Curfman McInnes     sequence as the usual matrix interface routines, since they
2059deebb3c3SLois Curfman McInnes     are intended to be accessed via the usual matrix interface
2060deebb3c3SLois Curfman McInnes     routines, e.g.,
2061a62d957aSLois Curfman McInnes $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
2062a62d957aSLois Curfman McInnes 
20634aa34b0aSBarry Smith     In particular each function MUST return an error code of 0 on success and
20644aa34b0aSBarry Smith     nonzero on failure.
20654aa34b0aSBarry Smith 
2066a62d957aSLois Curfman McInnes     Within each user-defined routine, the user should call
2067a62d957aSLois Curfman McInnes     MatShellGetContext() to obtain the user-defined context that was
2068a62d957aSLois Curfman McInnes     set by MatCreateShell().
2069a62d957aSLois Curfman McInnes 
2070b77ba244SStefano Zampini     Use MatSetOperation() to set an operation for any matrix type. For matrix product operations (i.e. MatMat, MatTransposeMat etc) use MatShellSetMatProductOperation()
2071b77ba244SStefano Zampini 
207295452b02SPatrick Sanan     Fortran Notes:
207395452b02SPatrick Sanan     For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
2074c4762a1bSJed Brown        generate a matrix. See src/mat/tests/ex120f.F
2075501d9185SBarry Smith 
2076b77ba244SStefano Zampini .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts(), MatShellSetMatProductOperation()
2077e51e0e81SBarry Smith @*/
2078789d8953SBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*g)(void))
2079e51e0e81SBarry Smith {
2080976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
2081273d9f13SBarry Smith 
20823a40ed3dSBarry Smith   PetscFunctionBegin;
20830700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2084ef5c7bd2SStefano Zampini   ierr = PetscTryMethod(mat,"MatShellSetOperation_C",(Mat,MatOperation,void (*)(void)),(mat,op,g));CHKERRQ(ierr);
20853a40ed3dSBarry Smith   PetscFunctionReturn(0);
2086e51e0e81SBarry Smith }
2087f0479e8cSBarry Smith 
2088d4bb536fSBarry Smith /*@C
2089d4bb536fSBarry Smith     MatShellGetOperation - Gets a matrix function for a shell matrix.
2090d4bb536fSBarry Smith 
2091c7fcc2eaSBarry Smith     Not Collective
2092c7fcc2eaSBarry Smith 
2093d4bb536fSBarry Smith     Input Parameters:
2094c7fcc2eaSBarry Smith +   mat - the shell matrix
2095c7fcc2eaSBarry Smith -   op - the name of the operation
2096d4bb536fSBarry Smith 
2097d4bb536fSBarry Smith     Output Parameter:
2098789d8953SBarry Smith .   g - the function that provides the operation.
2099d4bb536fSBarry Smith 
210015091d37SBarry Smith     Level: advanced
210115091d37SBarry Smith 
2102d4bb536fSBarry Smith     Notes:
2103e090d566SSatish Balay     See the file include/petscmat.h for a complete list of matrix
2104d4bb536fSBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
2105d4bb536fSBarry Smith     <OPERATION> is the name (in all capital letters) of the
2106d4bb536fSBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
2107d4bb536fSBarry Smith 
2108d4bb536fSBarry Smith     All user-provided functions have the same calling
2109d4bb536fSBarry Smith     sequence as the usual matrix interface routines, since they
2110d4bb536fSBarry Smith     are intended to be accessed via the usual matrix interface
2111d4bb536fSBarry Smith     routines, e.g.,
2112d4bb536fSBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
2113d4bb536fSBarry Smith 
2114d4bb536fSBarry Smith     Within each user-defined routine, the user should call
2115d4bb536fSBarry Smith     MatShellGetContext() to obtain the user-defined context that was
2116d4bb536fSBarry Smith     set by MatCreateShell().
2117d4bb536fSBarry Smith 
2118ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
2119d4bb536fSBarry Smith @*/
2120789d8953SBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**g)(void))
2121d4bb536fSBarry Smith {
2122976e8c5aSLisandro Dalcin   PetscErrorCode ierr;
2123273d9f13SBarry Smith 
21243a40ed3dSBarry Smith   PetscFunctionBegin;
21250700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2126789d8953SBarry Smith   ierr = PetscUseMethod(mat,"MatShellGetOperation_C",(Mat,MatOperation,void (**)(void)),(mat,op,g));CHKERRQ(ierr);
21273a40ed3dSBarry Smith   PetscFunctionReturn(0);
2128d4bb536fSBarry Smith }
2129b77ba244SStefano Zampini 
2130b77ba244SStefano Zampini /*@
2131b77ba244SStefano Zampini     MatIsShell - Inquires if a matrix is derived from MATSHELL
2132b77ba244SStefano Zampini 
2133b77ba244SStefano Zampini     Input Parameter:
2134b77ba244SStefano Zampini .   mat - the matrix
2135b77ba244SStefano Zampini 
2136b77ba244SStefano Zampini     Output Parameter:
2137b77ba244SStefano Zampini .   flg - the boolean value
2138b77ba244SStefano Zampini 
2139b77ba244SStefano Zampini     Level: developer
2140b77ba244SStefano Zampini 
2141b77ba244SStefano Zampini     Notes: in the future, we should allow the object type name to be changed still using the MatShell data structure for other matrices (i.e. MATTRANSPOSEMAT, MATSCHURCOMPLEMENT etc)
2142b77ba244SStefano Zampini 
2143b77ba244SStefano Zampini .seealso: MatCreateShell()
2144b77ba244SStefano Zampini @*/
2145b77ba244SStefano Zampini PetscErrorCode MatIsShell(Mat mat, PetscBool *flg)
2146b77ba244SStefano Zampini {
2147b77ba244SStefano Zampini   PetscFunctionBegin;
2148b77ba244SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2149b77ba244SStefano Zampini   PetscValidPointer(flg,2);
2150b77ba244SStefano Zampini   *flg = (PetscBool)(mat->ops->destroy == MatDestroy_Shell);
2151b77ba244SStefano Zampini   PetscFunctionReturn(0);
2152b77ba244SStefano Zampini }
2153