1be1d678aSKris Buschelman 2e51e0e81SBarry Smith /* 320563c6bSBarry Smith This provides a simple shell for Fortran (and C programmers) to 420563c6bSBarry Smith create a very simple matrix class for use with KSP without coding 5ed3cc1f0SBarry Smith much of anything. 6e51e0e81SBarry Smith */ 7e51e0e81SBarry Smith 8af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 945960306SStefano Zampini #include <petsc/private/vecimpl.h> 10e51e0e81SBarry Smith 1128f357e3SAlex Fikl struct _MatShellOps { 12976e8c5aSLisandro Dalcin /* 3 */ PetscErrorCode (*mult)(Mat,Vec,Vec); 13976e8c5aSLisandro Dalcin /* 5 */ PetscErrorCode (*multtranspose)(Mat,Vec,Vec); 14976e8c5aSLisandro Dalcin /* 17 */ PetscErrorCode (*getdiagonal)(Mat,Vec); 15976e8c5aSLisandro Dalcin /* 43 */ PetscErrorCode (*copy)(Mat,Mat,MatStructure); 16976e8c5aSLisandro Dalcin /* 60 */ PetscErrorCode (*destroy)(Mat); 1728f357e3SAlex Fikl }; 1828f357e3SAlex Fikl 1928f357e3SAlex Fikl typedef struct { 2028f357e3SAlex Fikl struct _MatShellOps ops[1]; 212205254eSKarl Rupp 22ef66eb69SBarry Smith PetscScalar vscale,vshift; 238fe8eb27SJed Brown Vec dshift; 248fe8eb27SJed Brown Vec left,right; 258fe8eb27SJed Brown Vec left_work,right_work; 265edf6516SJed Brown Vec left_add_work,right_add_work; 279f137db4SBarry Smith Mat axpy; 289f137db4SBarry Smith PetscScalar axpy_vscale; 290c0fd78eSBarry Smith PetscBool managescalingshifts; /* The user will manage the scaling and shifts for the MATSHELL, not the default */ 3045960306SStefano Zampini /* support for ZeroRows/Columns operations */ 3145960306SStefano Zampini IS zrows; 3245960306SStefano Zampini IS zcols; 3345960306SStefano Zampini Vec zvals; 3445960306SStefano Zampini Vec zvals_w; 3545960306SStefano Zampini VecScatter zvals_sct_r; 3645960306SStefano Zampini VecScatter zvals_sct_c; 3720563c6bSBarry Smith void *ctx; 3888cf3e7dSBarry Smith } Mat_Shell; 390c0fd78eSBarry Smith 4045960306SStefano Zampini 4145960306SStefano Zampini /* 4245960306SStefano Zampini Store and scale values on zeroed rows 4345960306SStefano Zampini xx = [x_1, 0], 0 on zeroed columns 4445960306SStefano Zampini */ 4545960306SStefano Zampini static PetscErrorCode MatShellPreZeroRight(Mat A,Vec x,Vec *xx) 4645960306SStefano Zampini { 4745960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 4845960306SStefano Zampini PetscErrorCode ierr; 4945960306SStefano Zampini 5045960306SStefano Zampini PetscFunctionBegin; 5145960306SStefano Zampini *xx = x; 5245960306SStefano Zampini if (shell->zrows) { 5345960306SStefano Zampini ierr = VecSet(shell->zvals_w,0.0);CHKERRQ(ierr); 5445960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5545960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5645960306SStefano Zampini ierr = VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);CHKERRQ(ierr); 5745960306SStefano Zampini } 5845960306SStefano Zampini if (shell->zcols) { 5945960306SStefano Zampini if (!shell->right_work) { 6045960306SStefano Zampini ierr = MatCreateVecs(A,&shell->right_work,NULL);CHKERRQ(ierr); 6145960306SStefano Zampini } 6245960306SStefano Zampini ierr = VecCopy(x,shell->right_work);CHKERRQ(ierr); 6345960306SStefano Zampini ierr = VecISSet(shell->right_work,shell->zcols,0.0);CHKERRQ(ierr); 6445960306SStefano Zampini *xx = shell->right_work; 6545960306SStefano Zampini } 6645960306SStefano Zampini PetscFunctionReturn(0); 6745960306SStefano Zampini } 6845960306SStefano Zampini 6945960306SStefano Zampini /* Insert properly diagonally scaled values stored in MatShellPreZeroRight */ 7045960306SStefano Zampini static PetscErrorCode MatShellPostZeroLeft(Mat A,Vec x) 7145960306SStefano Zampini { 7245960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 7345960306SStefano Zampini PetscErrorCode ierr; 7445960306SStefano Zampini 7545960306SStefano Zampini PetscFunctionBegin; 7645960306SStefano Zampini if (shell->zrows) { 7745960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7845960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7945960306SStefano Zampini } 8045960306SStefano Zampini PetscFunctionReturn(0); 8145960306SStefano Zampini } 8245960306SStefano Zampini 8345960306SStefano Zampini /* 8445960306SStefano Zampini Store and scale values on zeroed rows 8545960306SStefano Zampini xx = [x_1, 0], 0 on zeroed rows 8645960306SStefano Zampini */ 8745960306SStefano Zampini static PetscErrorCode MatShellPreZeroLeft(Mat A,Vec x,Vec *xx) 8845960306SStefano Zampini { 8945960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 9045960306SStefano Zampini PetscErrorCode ierr; 9145960306SStefano Zampini 9245960306SStefano Zampini PetscFunctionBegin; 9345960306SStefano Zampini *xx = NULL; 9445960306SStefano Zampini if (!shell->zrows) { 9545960306SStefano Zampini *xx = x; 9645960306SStefano Zampini } else { 9745960306SStefano Zampini if (!shell->left_work) { 9845960306SStefano Zampini ierr = MatCreateVecs(A,NULL,&shell->left_work);CHKERRQ(ierr); 9945960306SStefano Zampini } 10045960306SStefano Zampini ierr = VecCopy(x,shell->left_work);CHKERRQ(ierr); 10145960306SStefano Zampini ierr = VecSet(shell->zvals_w,0.0);CHKERRQ(ierr); 10245960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10345960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10445960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10545960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10645960306SStefano Zampini ierr = VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);CHKERRQ(ierr); 10745960306SStefano Zampini *xx = shell->left_work; 10845960306SStefano Zampini } 10945960306SStefano Zampini PetscFunctionReturn(0); 11045960306SStefano Zampini } 11145960306SStefano Zampini 11245960306SStefano Zampini /* Zero zero-columns contributions, sum contributions from properly scaled values stored in MatShellPreZeroLeft */ 11345960306SStefano Zampini static PetscErrorCode MatShellPostZeroRight(Mat A,Vec x) 11445960306SStefano Zampini { 11545960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 11645960306SStefano Zampini PetscErrorCode ierr; 11745960306SStefano Zampini 11845960306SStefano Zampini PetscFunctionBegin; 11945960306SStefano Zampini if (shell->zcols) { 12045960306SStefano Zampini ierr = VecISSet(x,shell->zcols,0.0);CHKERRQ(ierr); 12145960306SStefano Zampini } 12245960306SStefano Zampini if (shell->zrows) { 12345960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12445960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12545960306SStefano Zampini } 12645960306SStefano Zampini PetscFunctionReturn(0); 12745960306SStefano Zampini } 12845960306SStefano Zampini 1298fe8eb27SJed Brown /* 1300c0fd78eSBarry Smith xx = diag(left)*x 1318fe8eb27SJed Brown */ 1328fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 1338fe8eb27SJed Brown { 1348fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1358fe8eb27SJed Brown PetscErrorCode ierr; 1368fe8eb27SJed Brown 1378fe8eb27SJed Brown PetscFunctionBegin; 1380298fd71SBarry Smith *xx = NULL; 1398fe8eb27SJed Brown if (!shell->left) { 1408fe8eb27SJed Brown *xx = x; 1418fe8eb27SJed Brown } else { 1428fe8eb27SJed Brown if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);} 1438fe8eb27SJed Brown ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr); 1448fe8eb27SJed Brown *xx = shell->left_work; 1458fe8eb27SJed Brown } 1468fe8eb27SJed Brown PetscFunctionReturn(0); 1478fe8eb27SJed Brown } 1488fe8eb27SJed Brown 1490c0fd78eSBarry Smith /* 1500c0fd78eSBarry Smith xx = diag(right)*x 1510c0fd78eSBarry Smith */ 1528fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 1538fe8eb27SJed Brown { 1548fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1558fe8eb27SJed Brown PetscErrorCode ierr; 1568fe8eb27SJed Brown 1578fe8eb27SJed Brown PetscFunctionBegin; 1580298fd71SBarry Smith *xx = NULL; 1598fe8eb27SJed Brown if (!shell->right) { 1608fe8eb27SJed Brown *xx = x; 1618fe8eb27SJed Brown } else { 1628fe8eb27SJed Brown if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);} 1638fe8eb27SJed Brown ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr); 1648fe8eb27SJed Brown *xx = shell->right_work; 1658fe8eb27SJed Brown } 1668fe8eb27SJed Brown PetscFunctionReturn(0); 1678fe8eb27SJed Brown } 1688fe8eb27SJed Brown 1690c0fd78eSBarry Smith /* 1700c0fd78eSBarry Smith x = diag(left)*x 1710c0fd78eSBarry Smith */ 1728fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 1738fe8eb27SJed Brown { 1748fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1758fe8eb27SJed Brown PetscErrorCode ierr; 1768fe8eb27SJed Brown 1778fe8eb27SJed Brown PetscFunctionBegin; 1788fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);} 1798fe8eb27SJed Brown PetscFunctionReturn(0); 1808fe8eb27SJed Brown } 1818fe8eb27SJed Brown 1820c0fd78eSBarry Smith /* 1830c0fd78eSBarry Smith x = diag(right)*x 1840c0fd78eSBarry Smith */ 1858fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 1868fe8eb27SJed Brown { 1878fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1888fe8eb27SJed Brown PetscErrorCode ierr; 1898fe8eb27SJed Brown 1908fe8eb27SJed Brown PetscFunctionBegin; 1918fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);} 1928fe8eb27SJed Brown PetscFunctionReturn(0); 1938fe8eb27SJed Brown } 1948fe8eb27SJed Brown 1950c0fd78eSBarry Smith /* 1960c0fd78eSBarry Smith Y = vscale*Y + diag(dshift)*X + vshift*X 1970c0fd78eSBarry Smith 1980c0fd78eSBarry Smith On input Y already contains A*x 1990c0fd78eSBarry Smith */ 2008fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 2018fe8eb27SJed Brown { 2028fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2038fe8eb27SJed Brown PetscErrorCode ierr; 2048fe8eb27SJed Brown 2058fe8eb27SJed Brown PetscFunctionBegin; 2068fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 2078fe8eb27SJed Brown PetscInt i,m; 2088fe8eb27SJed Brown const PetscScalar *x,*d; 2098fe8eb27SJed Brown PetscScalar *y; 2108fe8eb27SJed Brown ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr); 2118fe8eb27SJed Brown ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr); 2128fe8eb27SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 2138fe8eb27SJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 2148fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 2158fe8eb27SJed Brown ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr); 2168fe8eb27SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 2178fe8eb27SJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 2180c0fd78eSBarry Smith } else { 219027c5fb5SJed Brown ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr); 2208fe8eb27SJed Brown } 221d4c7de66SBarry Smith if (shell->vshift != 0.0) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */ 2228fe8eb27SJed Brown PetscFunctionReturn(0); 2238fe8eb27SJed Brown } 224e51e0e81SBarry Smith 225789d8953SBarry Smith PetscErrorCode MatShellGetContext_Shell(Mat mat,void *ctx) 226789d8953SBarry Smith { 227789d8953SBarry Smith PetscFunctionBegin; 228789d8953SBarry Smith *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 229789d8953SBarry Smith PetscFunctionReturn(0); 230789d8953SBarry Smith } 231789d8953SBarry Smith 2329d225801SJed Brown /*@ 233a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 234b4fd4287SBarry Smith 235c7fcc2eaSBarry Smith Not Collective 236c7fcc2eaSBarry Smith 237b4fd4287SBarry Smith Input Parameter: 238b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 239b4fd4287SBarry Smith 240b4fd4287SBarry Smith Output Parameter: 241b4fd4287SBarry Smith . ctx - the user provided context 242b4fd4287SBarry Smith 24315091d37SBarry Smith Level: advanced 24415091d37SBarry Smith 24595452b02SPatrick Sanan Fortran Notes: 24695452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 247daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 248a62d957aSLois Curfman McInnes 249ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 2500bc0a719SBarry Smith @*/ 2518fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 252b4fd4287SBarry Smith { 253dfbe8321SBarry Smith PetscErrorCode ierr; 254273d9f13SBarry Smith 2553a40ed3dSBarry Smith PetscFunctionBegin; 2560700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2574482741eSBarry Smith PetscValidPointer(ctx,2); 258789d8953SBarry Smith ierr = PetscUseMethod(mat,"MatShellGetContext_C",(Mat,void*),(mat,ctx));CHKERRQ(ierr); 2593a40ed3dSBarry Smith PetscFunctionReturn(0); 260b4fd4287SBarry Smith } 261b4fd4287SBarry Smith 26245960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Local_Shell(Mat mat,PetscInt nr,PetscInt rows[],PetscInt nc,PetscInt cols[],PetscScalar diag,PetscBool rc) 26345960306SStefano Zampini { 26445960306SStefano Zampini PetscErrorCode ierr; 26545960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)mat->data; 26645960306SStefano Zampini Vec x = NULL,b = NULL; 26745960306SStefano Zampini IS is1, is2; 26845960306SStefano Zampini const PetscInt *ridxs; 26945960306SStefano Zampini PetscInt *idxs,*gidxs; 27045960306SStefano Zampini PetscInt cum,rst,cst,i; 27145960306SStefano Zampini 27245960306SStefano Zampini PetscFunctionBegin; 27345960306SStefano Zampini if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation"); 27445960306SStefano Zampini if (!shell->zvals) { 27545960306SStefano Zampini ierr = MatCreateVecs(mat,NULL,&shell->zvals);CHKERRQ(ierr); 27645960306SStefano Zampini } 27745960306SStefano Zampini if (!shell->zvals_w) { 27845960306SStefano Zampini ierr = VecDuplicate(shell->zvals,&shell->zvals_w);CHKERRQ(ierr); 27945960306SStefano Zampini } 28045960306SStefano Zampini ierr = MatGetOwnershipRange(mat,&rst,NULL);CHKERRQ(ierr); 28145960306SStefano Zampini ierr = MatGetOwnershipRangeColumn(mat,&cst,NULL);CHKERRQ(ierr); 28245960306SStefano Zampini 28345960306SStefano Zampini /* Expand/create index set of zeroed rows */ 28445960306SStefano Zampini ierr = PetscMalloc1(nr,&idxs);CHKERRQ(ierr); 28545960306SStefano Zampini for (i = 0; i < nr; i++) idxs[i] = rows[i] + rst; 28645960306SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,nr,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr); 28745960306SStefano Zampini ierr = ISSort(is1);CHKERRQ(ierr); 28845960306SStefano Zampini ierr = VecISSet(shell->zvals,is1,diag);CHKERRQ(ierr); 28945960306SStefano Zampini if (shell->zrows) { 29045960306SStefano Zampini ierr = ISSum(shell->zrows,is1,&is2);CHKERRQ(ierr); 29145960306SStefano Zampini ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr); 29245960306SStefano Zampini ierr = ISDestroy(&is1);CHKERRQ(ierr); 29345960306SStefano Zampini shell->zrows = is2; 29445960306SStefano Zampini } else shell->zrows = is1; 29545960306SStefano Zampini 29645960306SStefano Zampini /* Create scatters for diagonal values communications */ 29745960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr); 29845960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr); 29945960306SStefano Zampini 30045960306SStefano Zampini /* row scatter: from/to left vector */ 30145960306SStefano Zampini ierr = MatCreateVecs(mat,&x,&b);CHKERRQ(ierr); 30245960306SStefano Zampini ierr = VecScatterCreate(b,shell->zrows,shell->zvals_w,shell->zrows,&shell->zvals_sct_r);CHKERRQ(ierr); 30345960306SStefano Zampini 30445960306SStefano Zampini /* col scatter: from right vector to left vector */ 30545960306SStefano Zampini ierr = ISGetIndices(shell->zrows,&ridxs);CHKERRQ(ierr); 30645960306SStefano Zampini ierr = ISGetLocalSize(shell->zrows,&nr);CHKERRQ(ierr); 30745960306SStefano Zampini ierr = PetscMalloc1(nr,&gidxs);CHKERRQ(ierr); 30845960306SStefano Zampini for (i = 0, cum = 0; i < nr; i++) { 30945960306SStefano Zampini if (ridxs[i] >= mat->cmap->N) continue; 31045960306SStefano Zampini gidxs[cum] = ridxs[i]; 31145960306SStefano Zampini cum++; 31245960306SStefano Zampini } 31345960306SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)mat),cum,gidxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr); 31445960306SStefano Zampini ierr = VecScatterCreate(x,is1,shell->zvals_w,is1,&shell->zvals_sct_c);CHKERRQ(ierr); 31545960306SStefano Zampini ierr = ISDestroy(&is1);CHKERRQ(ierr); 31645960306SStefano Zampini ierr = VecDestroy(&x);CHKERRQ(ierr); 31745960306SStefano Zampini ierr = VecDestroy(&b);CHKERRQ(ierr); 31845960306SStefano Zampini 31945960306SStefano Zampini /* Expand/create index set of zeroed columns */ 32045960306SStefano Zampini if (rc) { 32145960306SStefano Zampini ierr = PetscMalloc1(nc,&idxs);CHKERRQ(ierr); 32245960306SStefano Zampini for (i = 0; i < nc; i++) idxs[i] = cols[i] + cst; 32345960306SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,nc,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr); 32445960306SStefano Zampini ierr = ISSort(is1);CHKERRQ(ierr); 32545960306SStefano Zampini if (shell->zcols) { 32645960306SStefano Zampini ierr = ISSum(shell->zcols,is1,&is2);CHKERRQ(ierr); 32745960306SStefano Zampini ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr); 32845960306SStefano Zampini ierr = ISDestroy(&is1);CHKERRQ(ierr); 32945960306SStefano Zampini shell->zcols = is2; 33045960306SStefano Zampini } else shell->zcols = is1; 33145960306SStefano Zampini } 33245960306SStefano Zampini PetscFunctionReturn(0); 33345960306SStefano Zampini } 33445960306SStefano Zampini 33545960306SStefano Zampini static PetscErrorCode MatZeroRows_Shell(Mat mat,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 33645960306SStefano Zampini { 33745960306SStefano Zampini PetscInt nr, *lrows; 33845960306SStefano Zampini PetscErrorCode ierr; 33945960306SStefano Zampini 34045960306SStefano Zampini PetscFunctionBegin; 34145960306SStefano Zampini if (x && b) { 34245960306SStefano Zampini Vec xt; 34345960306SStefano Zampini PetscScalar *vals; 34445960306SStefano Zampini PetscInt *gcols,i,st,nl,nc; 34545960306SStefano Zampini 34645960306SStefano Zampini ierr = PetscMalloc1(n,&gcols);CHKERRQ(ierr); 34745960306SStefano Zampini for (i = 0, nc = 0; i < n; i++) if (rows[i] < mat->cmap->N) gcols[nc++] = rows[i]; 34845960306SStefano Zampini 34945960306SStefano Zampini ierr = MatCreateVecs(mat,&xt,NULL);CHKERRQ(ierr); 35045960306SStefano Zampini ierr = VecCopy(x,xt);CHKERRQ(ierr); 35145960306SStefano Zampini ierr = PetscCalloc1(nc,&vals);CHKERRQ(ierr); 35245960306SStefano Zampini ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */ 35345960306SStefano Zampini ierr = PetscFree(vals);CHKERRQ(ierr); 35445960306SStefano Zampini ierr = VecAssemblyBegin(xt);CHKERRQ(ierr); 35545960306SStefano Zampini ierr = VecAssemblyEnd(xt);CHKERRQ(ierr); 35645960306SStefano Zampini ierr = VecAYPX(xt,-1.0,x);CHKERRQ(ierr); /* xt = [0, x2] */ 35745960306SStefano Zampini 35845960306SStefano Zampini ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr); 35945960306SStefano Zampini ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr); 36045960306SStefano Zampini ierr = VecGetArray(xt,&vals);CHKERRQ(ierr); 36145960306SStefano Zampini for (i = 0; i < nl; i++) { 36245960306SStefano Zampini PetscInt g = i + st; 36345960306SStefano Zampini if (g > mat->rmap->N) continue; 36445960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 36545960306SStefano Zampini ierr = VecSetValue(b,g,diag*vals[i],INSERT_VALUES);CHKERRQ(ierr); 36645960306SStefano Zampini } 36745960306SStefano Zampini ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr); 36845960306SStefano Zampini ierr = VecAssemblyBegin(b);CHKERRQ(ierr); 36945960306SStefano Zampini ierr = VecAssemblyEnd(b);CHKERRQ(ierr); /* b = [b1, x2 * diag] */ 37045960306SStefano Zampini ierr = VecDestroy(&xt);CHKERRQ(ierr); 37145960306SStefano Zampini ierr = PetscFree(gcols);CHKERRQ(ierr); 37245960306SStefano Zampini } 37345960306SStefano Zampini ierr = PetscLayoutMapLocal(mat->rmap,n,rows,&nr,&lrows,NULL);CHKERRQ(ierr); 37445960306SStefano Zampini ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,0,NULL,diag,PETSC_FALSE);CHKERRQ(ierr); 37545960306SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 37645960306SStefano Zampini PetscFunctionReturn(0); 37745960306SStefano Zampini } 37845960306SStefano Zampini 37945960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Shell(Mat mat,PetscInt n,const PetscInt rowscols[],PetscScalar diag,Vec x,Vec b) 38045960306SStefano Zampini { 38145960306SStefano Zampini PetscInt *lrows, *lcols; 38245960306SStefano Zampini PetscInt nr, nc; 38345960306SStefano Zampini PetscBool congruent; 38445960306SStefano Zampini PetscErrorCode ierr; 38545960306SStefano Zampini 38645960306SStefano Zampini PetscFunctionBegin; 38745960306SStefano Zampini if (x && b) { 38845960306SStefano Zampini Vec xt, bt; 38945960306SStefano Zampini PetscScalar *vals; 39045960306SStefano Zampini PetscInt *grows,*gcols,i,st,nl; 39145960306SStefano Zampini 39245960306SStefano Zampini ierr = PetscMalloc2(n,&grows,n,&gcols);CHKERRQ(ierr); 39345960306SStefano Zampini for (i = 0, nr = 0; i < n; i++) if (rowscols[i] < mat->rmap->N) grows[nr++] = rowscols[i]; 39445960306SStefano Zampini for (i = 0, nc = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) gcols[nc++] = rowscols[i]; 39545960306SStefano Zampini ierr = PetscCalloc1(n,&vals);CHKERRQ(ierr); 39645960306SStefano Zampini 39745960306SStefano Zampini ierr = MatCreateVecs(mat,&xt,&bt);CHKERRQ(ierr); 39845960306SStefano Zampini ierr = VecCopy(x,xt);CHKERRQ(ierr); 39945960306SStefano Zampini ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */ 40045960306SStefano Zampini ierr = VecAssemblyBegin(xt);CHKERRQ(ierr); 40145960306SStefano Zampini ierr = VecAssemblyEnd(xt);CHKERRQ(ierr); 40245960306SStefano Zampini ierr = VecAXPY(xt,-1.0,x);CHKERRQ(ierr); /* xt = [0, -x2] */ 40345960306SStefano Zampini ierr = MatMult(mat,xt,bt);CHKERRQ(ierr); /* bt = [-A12*x2,-A22*x2] */ 40445960306SStefano Zampini ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* bt = [-A12*x2,0] */ 40545960306SStefano Zampini ierr = VecAssemblyBegin(bt);CHKERRQ(ierr); 40645960306SStefano Zampini ierr = VecAssemblyEnd(bt);CHKERRQ(ierr); 40745960306SStefano Zampini ierr = VecAXPY(b,1.0,bt);CHKERRQ(ierr); /* b = [b1 - A12*x2, b2] */ 40845960306SStefano Zampini ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* b = [b1 - A12*x2, 0] */ 40945960306SStefano Zampini ierr = VecAssemblyBegin(bt);CHKERRQ(ierr); 41045960306SStefano Zampini ierr = VecAssemblyEnd(bt);CHKERRQ(ierr); 41145960306SStefano Zampini ierr = PetscFree(vals);CHKERRQ(ierr); 41245960306SStefano Zampini 41345960306SStefano Zampini ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr); 41445960306SStefano Zampini ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr); 41545960306SStefano Zampini ierr = VecGetArray(xt,&vals);CHKERRQ(ierr); 41645960306SStefano Zampini for (i = 0; i < nl; i++) { 41745960306SStefano Zampini PetscInt g = i + st; 41845960306SStefano Zampini if (g > mat->rmap->N) continue; 41945960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 42045960306SStefano Zampini ierr = VecSetValue(b,g,-diag*vals[i],INSERT_VALUES);CHKERRQ(ierr); 42145960306SStefano Zampini } 42245960306SStefano Zampini ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr); 42345960306SStefano Zampini ierr = VecAssemblyBegin(b);CHKERRQ(ierr); 42445960306SStefano Zampini ierr = VecAssemblyEnd(b);CHKERRQ(ierr); /* b = [b1 - A12*x2, x2 * diag] */ 42545960306SStefano Zampini ierr = VecDestroy(&xt);CHKERRQ(ierr); 42645960306SStefano Zampini ierr = VecDestroy(&bt);CHKERRQ(ierr); 42745960306SStefano Zampini ierr = PetscFree2(grows,gcols);CHKERRQ(ierr); 42845960306SStefano Zampini } 42945960306SStefano Zampini ierr = PetscLayoutMapLocal(mat->rmap,n,rowscols,&nr,&lrows,NULL);CHKERRQ(ierr); 43045960306SStefano Zampini ierr = MatHasCongruentLayouts(mat,&congruent);CHKERRQ(ierr); 43145960306SStefano Zampini if (congruent) { 43245960306SStefano Zampini nc = nr; 43345960306SStefano Zampini lcols = lrows; 43445960306SStefano Zampini } else { /* MatZeroRowsColumns implicitly assumes the rowscols indices are for a square matrix, here we handle a more general case */ 43545960306SStefano Zampini PetscInt i,nt,*t; 43645960306SStefano Zampini 43745960306SStefano Zampini ierr = PetscMalloc1(n,&t);CHKERRQ(ierr); 43845960306SStefano Zampini for (i = 0, nt = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) t[nt++] = rowscols[i]; 43945960306SStefano Zampini ierr = PetscLayoutMapLocal(mat->cmap,nt,t,&nc,&lcols,NULL);CHKERRQ(ierr); 44045960306SStefano Zampini ierr = PetscFree(t);CHKERRQ(ierr); 44145960306SStefano Zampini } 44245960306SStefano Zampini ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,nc,lcols,diag,PETSC_TRUE);CHKERRQ(ierr); 44345960306SStefano Zampini if (!congruent) { 44445960306SStefano Zampini ierr = PetscFree(lcols);CHKERRQ(ierr); 44545960306SStefano Zampini } 44645960306SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 44745960306SStefano Zampini PetscFunctionReturn(0); 44845960306SStefano Zampini } 44945960306SStefano Zampini 450dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 451e51e0e81SBarry Smith { 452dfbe8321SBarry Smith PetscErrorCode ierr; 453bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 454ed3cc1f0SBarry Smith 4553a40ed3dSBarry Smith PetscFunctionBegin; 45628f357e3SAlex Fikl if (shell->ops->destroy) { 45728f357e3SAlex Fikl ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr); 458bf0cc555SLisandro Dalcin } 4590c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 4600c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 4610c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 4628fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 4638fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 4645edf6516SJed Brown ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr); 4655edf6516SJed Brown ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr); 4669f137db4SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 46745960306SStefano Zampini ierr = VecDestroy(&shell->zvals_w);CHKERRQ(ierr); 46845960306SStefano Zampini ierr = VecDestroy(&shell->zvals);CHKERRQ(ierr); 46945960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr); 47045960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr); 47145960306SStefano Zampini ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr); 47245960306SStefano Zampini ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr); 473bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 474789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellGetContext_C",NULL);CHKERRQ(ierr); 475789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetContext_C",NULL);CHKERRQ(ierr); 476*db77db73SJeremy L Thompson ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetVecType_C",NULL);CHKERRQ(ierr); 477789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetManageScalingShifts_C",NULL);CHKERRQ(ierr); 478789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetOperation_C",NULL);CHKERRQ(ierr); 479789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellGetOperation_C",NULL);CHKERRQ(ierr); 4803a40ed3dSBarry Smith PetscFunctionReturn(0); 481e51e0e81SBarry Smith } 482e51e0e81SBarry Smith 4837fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str) 4847fabda3fSAlex Fikl { 48528f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data; 4867fabda3fSAlex Fikl PetscErrorCode ierr; 487cb8c8a77SBarry Smith PetscBool matflg; 4887fabda3fSAlex Fikl 4897fabda3fSAlex Fikl PetscFunctionBegin; 49028f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr); 491cb8c8a77SBarry Smith if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL"); 4927fabda3fSAlex Fikl 493cb8c8a77SBarry Smith ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 494cb8c8a77SBarry Smith ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr); 4957fabda3fSAlex Fikl 496cb8c8a77SBarry Smith if (shellA->ops->copy) { 49728f357e3SAlex Fikl ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr); 498cb8c8a77SBarry Smith } 4997fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 5007fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 5010c0fd78eSBarry Smith if (shellA->dshift) { 5020c0fd78eSBarry Smith if (!shellB->dshift) { 5030c0fd78eSBarry Smith ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr); 5047fabda3fSAlex Fikl } 5050c0fd78eSBarry Smith ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr); 5067fabda3fSAlex Fikl } else { 5079f137db4SBarry Smith ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr); 5087fabda3fSAlex Fikl } 5090c0fd78eSBarry Smith if (shellA->left) { 5100c0fd78eSBarry Smith if (!shellB->left) { 5110c0fd78eSBarry Smith ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr); 5127fabda3fSAlex Fikl } 5130c0fd78eSBarry Smith ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr); 5147fabda3fSAlex Fikl } else { 5159f137db4SBarry Smith ierr = VecDestroy(&shellB->left);CHKERRQ(ierr); 5167fabda3fSAlex Fikl } 5170c0fd78eSBarry Smith if (shellA->right) { 5180c0fd78eSBarry Smith if (!shellB->right) { 5190c0fd78eSBarry Smith ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr); 5207fabda3fSAlex Fikl } 5210c0fd78eSBarry Smith ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr); 5227fabda3fSAlex Fikl } else { 5239f137db4SBarry Smith ierr = VecDestroy(&shellB->right);CHKERRQ(ierr); 5247fabda3fSAlex Fikl } 52540e381c3SBarry Smith ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr); 52640e381c3SBarry Smith if (shellA->axpy) { 52740e381c3SBarry Smith ierr = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr); 52840e381c3SBarry Smith shellB->axpy = shellA->axpy; 52940e381c3SBarry Smith shellB->axpy_vscale = shellA->axpy_vscale; 53040e381c3SBarry Smith } 53145960306SStefano Zampini if (shellA->zrows) { 53245960306SStefano Zampini ierr = ISDuplicate(shellA->zrows,&shellB->zrows);CHKERRQ(ierr); 53345960306SStefano Zampini if (shellA->zcols) { 53445960306SStefano Zampini ierr = ISDuplicate(shellA->zcols,&shellB->zcols);CHKERRQ(ierr); 53545960306SStefano Zampini } 53645960306SStefano Zampini ierr = VecDuplicate(shellA->zvals,&shellB->zvals);CHKERRQ(ierr); 53745960306SStefano Zampini ierr = VecCopy(shellA->zvals,shellB->zvals);CHKERRQ(ierr); 53845960306SStefano Zampini ierr = VecDuplicate(shellA->zvals_w,&shellB->zvals_w);CHKERRQ(ierr); 53945960306SStefano Zampini ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_r);CHKERRQ(ierr); 54045960306SStefano Zampini ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_c);CHKERRQ(ierr); 54145960306SStefano Zampini shellB->zvals_sct_r = shellA->zvals_sct_r; 54245960306SStefano Zampini shellB->zvals_sct_c = shellA->zvals_sct_c; 54345960306SStefano Zampini } 5447fabda3fSAlex Fikl PetscFunctionReturn(0); 5457fabda3fSAlex Fikl } 5467fabda3fSAlex Fikl 547cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M) 548cb8c8a77SBarry Smith { 549cb8c8a77SBarry Smith PetscErrorCode ierr; 550cb8c8a77SBarry Smith void *ctx; 551cb8c8a77SBarry Smith 552cb8c8a77SBarry Smith PetscFunctionBegin; 553cb8c8a77SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 554cb8c8a77SBarry Smith ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr); 555a4b1380bSStefano Zampini if (op != MAT_DO_NOT_COPY_VALUES) { 556cb8c8a77SBarry Smith ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 557a4b1380bSStefano Zampini } 558cb8c8a77SBarry Smith PetscFunctionReturn(0); 559cb8c8a77SBarry Smith } 560cb8c8a77SBarry Smith 561dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 562ef66eb69SBarry Smith { 563ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 564dfbe8321SBarry Smith PetscErrorCode ierr; 56525578ef6SJed Brown Vec xx; 566e3f487b0SBarry Smith PetscObjectState instate,outstate; 567ef66eb69SBarry Smith 568ef66eb69SBarry Smith PetscFunctionBegin; 569976e8c5aSLisandro Dalcin if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL"); 57045960306SStefano Zampini ierr = MatShellPreZeroRight(A,x,&xx);CHKERRQ(ierr); 57145960306SStefano Zampini ierr = MatShellPreScaleRight(A,xx,&xx);CHKERRQ(ierr); 57245960306SStefano Zampini ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 57328f357e3SAlex Fikl ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr); 574e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 575e3f487b0SBarry Smith if (instate == outstate) { 576e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 577e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 578e3f487b0SBarry Smith } 5798fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 5808fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 58145960306SStefano Zampini ierr = MatShellPostZeroLeft(A,y);CHKERRQ(ierr); 5829f137db4SBarry Smith 5839f137db4SBarry Smith if (shell->axpy) { 5849f137db4SBarry Smith if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);} 5859f137db4SBarry Smith ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr); 5869f137db4SBarry Smith ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 5879f137db4SBarry Smith } 588ef66eb69SBarry Smith PetscFunctionReturn(0); 589ef66eb69SBarry Smith } 590ef66eb69SBarry Smith 5915edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 5925edf6516SJed Brown { 5935edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 5945edf6516SJed Brown PetscErrorCode ierr; 5955edf6516SJed Brown 5965edf6516SJed Brown PetscFunctionBegin; 5975edf6516SJed Brown if (y == z) { 5985edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 5995edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 600b8430d49SBarry Smith ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr); 6015edf6516SJed Brown } else { 6025edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 6035edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 6045edf6516SJed Brown } 6055edf6516SJed Brown PetscFunctionReturn(0); 6065edf6516SJed Brown } 6075edf6516SJed Brown 60818be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 60918be62a5SSatish Balay { 61018be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 61118be62a5SSatish Balay PetscErrorCode ierr; 61225578ef6SJed Brown Vec xx; 613e3f487b0SBarry Smith PetscObjectState instate,outstate; 61418be62a5SSatish Balay 61518be62a5SSatish Balay PetscFunctionBegin; 61645960306SStefano Zampini ierr = MatShellPreZeroLeft(A,x,&xx);CHKERRQ(ierr); 61745960306SStefano Zampini ierr = MatShellPreScaleLeft(A,xx,&xx);CHKERRQ(ierr); 618e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 6190c0fd78eSBarry Smith if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL"); 62028f357e3SAlex Fikl ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr); 621e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 622e3f487b0SBarry Smith if (instate == outstate) { 623e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 624e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 625e3f487b0SBarry Smith } 6268fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 6278fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 62845960306SStefano Zampini ierr = MatShellPostZeroRight(A,y);CHKERRQ(ierr); 629350ec83bSStefano Zampini 630350ec83bSStefano Zampini if (shell->axpy) { 631350ec83bSStefano Zampini if (!shell->right_work) {ierr = MatCreateVecs(A,NULL,&shell->right_work);CHKERRQ(ierr);} 632350ec83bSStefano Zampini ierr = MatMultTranspose(shell->axpy,x,shell->right_work);CHKERRQ(ierr); 633350ec83bSStefano Zampini ierr = VecAXPY(y,shell->axpy_vscale,shell->right_work);CHKERRQ(ierr); 634350ec83bSStefano Zampini } 63518be62a5SSatish Balay PetscFunctionReturn(0); 63618be62a5SSatish Balay } 63718be62a5SSatish Balay 6385edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 6395edf6516SJed Brown { 6405edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 6415edf6516SJed Brown PetscErrorCode ierr; 6425edf6516SJed Brown 6435edf6516SJed Brown PetscFunctionBegin; 6445edf6516SJed Brown if (y == z) { 6455edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 6465edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 647dac989eaS“Marek ierr = VecAXPY(z,1.0,shell->left_add_work);CHKERRQ(ierr); 6485edf6516SJed Brown } else { 6495edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 6505edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 6515edf6516SJed Brown } 6525edf6516SJed Brown PetscFunctionReturn(0); 6535edf6516SJed Brown } 6545edf6516SJed Brown 6550c0fd78eSBarry Smith /* 6560c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 6570c0fd78eSBarry Smith */ 65818be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 65918be62a5SSatish Balay { 66018be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 66118be62a5SSatish Balay PetscErrorCode ierr; 66218be62a5SSatish Balay 66318be62a5SSatish Balay PetscFunctionBegin; 6640c0fd78eSBarry Smith if (shell->ops->getdiagonal) { 66528f357e3SAlex Fikl ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr); 666305211d5SBarry 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,...)"); 66718be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 6688fe8eb27SJed Brown if (shell->dshift) { 6690c0fd78eSBarry Smith ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr); 67018be62a5SSatish Balay } 6710c0fd78eSBarry Smith ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 6728fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 6738fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 67445960306SStefano Zampini if (shell->zrows) { 67545960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 67645960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 67745960306SStefano Zampini } 67881c02519SBarry Smith if (shell->axpy) { 67981c02519SBarry Smith if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);} 68081c02519SBarry Smith ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr); 68181c02519SBarry Smith ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 68281c02519SBarry Smith } 68318be62a5SSatish Balay PetscFunctionReturn(0); 68418be62a5SSatish Balay } 68518be62a5SSatish Balay 686f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 687ef66eb69SBarry Smith { 688ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 6898fe8eb27SJed Brown PetscErrorCode ierr; 690789d8953SBarry Smith PetscBool flg; 691b24ad042SBarry Smith 692ef66eb69SBarry Smith PetscFunctionBegin; 693789d8953SBarry Smith ierr = MatHasCongruentLayouts(Y,&flg);CHKERRQ(ierr); 694789d8953SBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot shift shell matrix if it is not congruent"); 6950c0fd78eSBarry Smith if (shell->left || shell->right) { 6968fe8eb27SJed Brown if (!shell->dshift) { 6970c0fd78eSBarry Smith ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr); 6980c0fd78eSBarry Smith ierr = VecSet(shell->dshift,a);CHKERRQ(ierr); 6990c0fd78eSBarry Smith } else { 7000c0fd78eSBarry Smith if (shell->left) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 7010c0fd78eSBarry Smith if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 7029f137db4SBarry Smith ierr = VecShift(shell->dshift,a);CHKERRQ(ierr); 7030c0fd78eSBarry Smith } 7048fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 7058fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 7068fe8eb27SJed Brown } else shell->vshift += a; 70745960306SStefano Zampini if (shell->zrows) { 70845960306SStefano Zampini ierr = VecShift(shell->zvals,a);CHKERRQ(ierr); 70945960306SStefano Zampini } 710ef66eb69SBarry Smith PetscFunctionReturn(0); 711ef66eb69SBarry Smith } 712ef66eb69SBarry Smith 713b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s) 7146464bf51SAlex Fikl { 7156464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 7166464bf51SAlex Fikl PetscErrorCode ierr; 7176464bf51SAlex Fikl 7186464bf51SAlex Fikl PetscFunctionBegin; 7190c0fd78eSBarry Smith if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);} 7200c0fd78eSBarry Smith if (shell->left || shell->right) { 72192fabd1fSBarry Smith if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);} 7229f137db4SBarry Smith if (shell->left && shell->right) { 7239f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 7249f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr); 7259f137db4SBarry Smith } else if (shell->left) { 7269f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 7279f137db4SBarry Smith } else { 7289f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr); 7299f137db4SBarry Smith } 730b253ae0bS“Marek ierr = VecAXPY(shell->dshift,s,shell->right_work);CHKERRQ(ierr); 7310c0fd78eSBarry Smith } else { 732b253ae0bS“Marek ierr = VecAXPY(shell->dshift,s,D);CHKERRQ(ierr); 733b253ae0bS“Marek } 734b253ae0bS“Marek PetscFunctionReturn(0); 735b253ae0bS“Marek } 736b253ae0bS“Marek 737b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 738b253ae0bS“Marek { 73945960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 740b253ae0bS“Marek Vec d; 741b253ae0bS“Marek PetscErrorCode ierr; 742789d8953SBarry Smith PetscBool flg; 743b253ae0bS“Marek 744b253ae0bS“Marek PetscFunctionBegin; 745789d8953SBarry Smith ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr); 746789d8953SBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot diagonal set or shift shell matrix if it is not congruent"); 747b253ae0bS“Marek if (ins == INSERT_VALUES) { 748b253ae0bS“Marek if (!A->ops->getdiagonal) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first"); 749b253ae0bS“Marek ierr = VecDuplicate(D,&d);CHKERRQ(ierr); 750b253ae0bS“Marek ierr = MatGetDiagonal(A,d);CHKERRQ(ierr); 751b253ae0bS“Marek ierr = MatDiagonalSet_Shell_Private(A,d,-1.);CHKERRQ(ierr); 752b253ae0bS“Marek ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr); 753b253ae0bS“Marek ierr = VecDestroy(&d);CHKERRQ(ierr); 75445960306SStefano Zampini if (shell->zrows) { 75545960306SStefano Zampini ierr = VecCopy(D,shell->zvals);CHKERRQ(ierr); 75645960306SStefano Zampini } 757b253ae0bS“Marek } else { 758b253ae0bS“Marek ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr); 75945960306SStefano Zampini if (shell->zrows) { 76045960306SStefano Zampini ierr = VecAXPY(shell->zvals,1.0,D);CHKERRQ(ierr); 76145960306SStefano Zampini } 7626464bf51SAlex Fikl } 7636464bf51SAlex Fikl PetscFunctionReturn(0); 7646464bf51SAlex Fikl } 7656464bf51SAlex Fikl 766f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 767ef66eb69SBarry Smith { 768ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 7698fe8eb27SJed Brown PetscErrorCode ierr; 770b24ad042SBarry Smith 771ef66eb69SBarry Smith PetscFunctionBegin; 772f4df32b1SMatthew Knepley shell->vscale *= a; 7730c0fd78eSBarry Smith shell->vshift *= a; 7742205254eSKarl Rupp if (shell->dshift) { 7752205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 7760c0fd78eSBarry Smith } 77781c02519SBarry Smith shell->axpy_vscale *= a; 77845960306SStefano Zampini if (shell->zrows) { 77945960306SStefano Zampini ierr = VecScale(shell->zvals,a);CHKERRQ(ierr); 78045960306SStefano Zampini } 7818fe8eb27SJed Brown PetscFunctionReturn(0); 78218be62a5SSatish Balay } 7838fe8eb27SJed Brown 7848fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 7858fe8eb27SJed Brown { 7868fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 7878fe8eb27SJed Brown PetscErrorCode ierr; 7888fe8eb27SJed Brown 7898fe8eb27SJed Brown PetscFunctionBegin; 79081c02519SBarry Smith if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation"); 7918fe8eb27SJed Brown if (left) { 7920c0fd78eSBarry Smith if (!shell->left) { 7930c0fd78eSBarry Smith ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr); 7948fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 7950c0fd78eSBarry Smith } else { 7960c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 79718be62a5SSatish Balay } 79845960306SStefano Zampini if (shell->zrows) { 79945960306SStefano Zampini ierr = VecPointwiseMult(shell->zvals,shell->zvals,left);CHKERRQ(ierr); 80045960306SStefano Zampini } 801ef66eb69SBarry Smith } 8028fe8eb27SJed Brown if (right) { 8030c0fd78eSBarry Smith if (!shell->right) { 8040c0fd78eSBarry Smith ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr); 8058fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 8060c0fd78eSBarry Smith } else { 8070c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 8088fe8eb27SJed Brown } 80945960306SStefano Zampini if (shell->zrows) { 81045960306SStefano Zampini if (!shell->left_work) { 81145960306SStefano Zampini ierr = MatCreateVecs(Y,NULL,&shell->left_work);CHKERRQ(ierr); 81245960306SStefano Zampini } 81345960306SStefano Zampini ierr = VecSet(shell->zvals_w,1.0);CHKERRQ(ierr); 81445960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 81545960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 81645960306SStefano Zampini ierr = VecPointwiseMult(shell->zvals,shell->zvals,shell->zvals_w);CHKERRQ(ierr); 81745960306SStefano Zampini } 8188fe8eb27SJed Brown } 819ef66eb69SBarry Smith PetscFunctionReturn(0); 820ef66eb69SBarry Smith } 821ef66eb69SBarry Smith 822dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 823ef66eb69SBarry Smith { 824ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 8250c0fd78eSBarry Smith PetscErrorCode ierr; 826ef66eb69SBarry Smith 827ef66eb69SBarry Smith PetscFunctionBegin; 8288fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 829ef66eb69SBarry Smith shell->vshift = 0.0; 830ef66eb69SBarry Smith shell->vscale = 1.0; 8310c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 8320c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 8330c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 83481c02519SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 83545960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr); 83645960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr); 83745960306SStefano Zampini ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr); 83845960306SStefano Zampini ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr); 839ef66eb69SBarry Smith } 840ef66eb69SBarry Smith PetscFunctionReturn(0); 841ef66eb69SBarry Smith } 842ef66eb69SBarry Smith 8433b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 8443b49f96aSBarry Smith { 8453b49f96aSBarry Smith PetscFunctionBegin; 8463b49f96aSBarry Smith *missing = PETSC_FALSE; 8473b49f96aSBarry Smith PetscFunctionReturn(0); 8483b49f96aSBarry Smith } 8493b49f96aSBarry Smith 8509f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str) 8519f137db4SBarry Smith { 8529f137db4SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 8539f137db4SBarry Smith PetscErrorCode ierr; 8549f137db4SBarry Smith 8559f137db4SBarry Smith PetscFunctionBegin; 8569f137db4SBarry Smith ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 8579f137db4SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 8589f137db4SBarry Smith shell->axpy = X; 8599f137db4SBarry Smith shell->axpy_vscale = a; 8609f137db4SBarry Smith PetscFunctionReturn(0); 8619f137db4SBarry Smith } 8629f137db4SBarry Smith 86309dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 86420563c6bSBarry Smith 0, 86520563c6bSBarry Smith 0, 8669f137db4SBarry Smith 0, 8670c0fd78eSBarry Smith /* 4*/ MatMultAdd_Shell, 8689f137db4SBarry Smith 0, 8690c0fd78eSBarry Smith MatMultTransposeAdd_Shell, 870b951964fSBarry Smith 0, 871b951964fSBarry Smith 0, 872b951964fSBarry Smith 0, 87397304618SKris Buschelman /*10*/ 0, 874b951964fSBarry Smith 0, 875b951964fSBarry Smith 0, 876b951964fSBarry Smith 0, 877b951964fSBarry Smith 0, 87897304618SKris Buschelman /*15*/ 0, 879b951964fSBarry Smith 0, 88000849d43SBarry Smith 0, 8818fe8eb27SJed Brown MatDiagonalScale_Shell, 882b951964fSBarry Smith 0, 88397304618SKris Buschelman /*20*/ 0, 884ef66eb69SBarry Smith MatAssemblyEnd_Shell, 885b951964fSBarry Smith 0, 886b951964fSBarry Smith 0, 88745960306SStefano Zampini /*24*/ MatZeroRows_Shell, 888b951964fSBarry Smith 0, 889b951964fSBarry Smith 0, 890b951964fSBarry Smith 0, 891b951964fSBarry Smith 0, 892d519adbfSMatthew Knepley /*29*/ 0, 893b951964fSBarry Smith 0, 894273d9f13SBarry Smith 0, 895b951964fSBarry Smith 0, 896b951964fSBarry Smith 0, 897cb8c8a77SBarry Smith /*34*/ MatDuplicate_Shell, 898b951964fSBarry Smith 0, 899b951964fSBarry Smith 0, 90009dc0095SBarry Smith 0, 90109dc0095SBarry Smith 0, 9029f137db4SBarry Smith /*39*/ MatAXPY_Shell, 90309dc0095SBarry Smith 0, 90409dc0095SBarry Smith 0, 90509dc0095SBarry Smith 0, 906cb8c8a77SBarry Smith MatCopy_Shell, 907d519adbfSMatthew Knepley /*44*/ 0, 908ef66eb69SBarry Smith MatScale_Shell, 909ef66eb69SBarry Smith MatShift_Shell, 9106464bf51SAlex Fikl MatDiagonalSet_Shell, 91145960306SStefano Zampini MatZeroRowsColumns_Shell, 912f73d5cc4SBarry Smith /*49*/ 0, 91309dc0095SBarry Smith 0, 91409dc0095SBarry Smith 0, 91509dc0095SBarry Smith 0, 91609dc0095SBarry Smith 0, 917d519adbfSMatthew Knepley /*54*/ 0, 91809dc0095SBarry Smith 0, 91909dc0095SBarry Smith 0, 92009dc0095SBarry Smith 0, 92109dc0095SBarry Smith 0, 922d519adbfSMatthew Knepley /*59*/ 0, 923b9b97703SBarry Smith MatDestroy_Shell, 92409dc0095SBarry Smith 0, 925251fad3fSStefano Zampini MatConvertFrom_Shell, 926273d9f13SBarry Smith 0, 927d519adbfSMatthew Knepley /*64*/ 0, 928273d9f13SBarry Smith 0, 929273d9f13SBarry Smith 0, 930273d9f13SBarry Smith 0, 931273d9f13SBarry Smith 0, 932d519adbfSMatthew Knepley /*69*/ 0, 933273d9f13SBarry Smith 0, 934c87e5d42SMatthew Knepley MatConvert_Shell, 935273d9f13SBarry Smith 0, 93697304618SKris Buschelman 0, 937d519adbfSMatthew Knepley /*74*/ 0, 93897304618SKris Buschelman 0, 93997304618SKris Buschelman 0, 94097304618SKris Buschelman 0, 94197304618SKris Buschelman 0, 942d519adbfSMatthew Knepley /*79*/ 0, 94397304618SKris Buschelman 0, 94497304618SKris Buschelman 0, 94597304618SKris Buschelman 0, 946865e5f61SKris Buschelman 0, 947d519adbfSMatthew Knepley /*84*/ 0, 948865e5f61SKris Buschelman 0, 949865e5f61SKris Buschelman 0, 950865e5f61SKris Buschelman 0, 951865e5f61SKris Buschelman 0, 952d519adbfSMatthew Knepley /*89*/ 0, 953865e5f61SKris Buschelman 0, 954865e5f61SKris Buschelman 0, 955865e5f61SKris Buschelman 0, 956865e5f61SKris Buschelman 0, 957d519adbfSMatthew Knepley /*94*/ 0, 958865e5f61SKris Buschelman 0, 959865e5f61SKris Buschelman 0, 9603964eb88SJed Brown 0, 9613964eb88SJed Brown 0, 9623964eb88SJed Brown /*99*/ 0, 9633964eb88SJed Brown 0, 9643964eb88SJed Brown 0, 9653964eb88SJed Brown 0, 9663964eb88SJed Brown 0, 9673964eb88SJed Brown /*104*/ 0, 9683964eb88SJed Brown 0, 9693964eb88SJed Brown 0, 9703964eb88SJed Brown 0, 9713964eb88SJed Brown 0, 9723964eb88SJed Brown /*109*/ 0, 9733964eb88SJed Brown 0, 9743964eb88SJed Brown 0, 9753964eb88SJed Brown 0, 9763b49f96aSBarry Smith MatMissingDiagonal_Shell, 9773964eb88SJed Brown /*114*/ 0, 9783964eb88SJed Brown 0, 9793964eb88SJed Brown 0, 9803964eb88SJed Brown 0, 9813964eb88SJed Brown 0, 9823964eb88SJed Brown /*119*/ 0, 9833964eb88SJed Brown 0, 9843964eb88SJed Brown 0, 9853964eb88SJed Brown 0, 9863964eb88SJed Brown 0, 9873964eb88SJed Brown /*124*/ 0, 9883964eb88SJed Brown 0, 9893964eb88SJed Brown 0, 9903964eb88SJed Brown 0, 9913964eb88SJed Brown 0, 9923964eb88SJed Brown /*129*/ 0, 9933964eb88SJed Brown 0, 9943964eb88SJed Brown 0, 9953964eb88SJed Brown 0, 9963964eb88SJed Brown 0, 9973964eb88SJed Brown /*134*/ 0, 9983964eb88SJed Brown 0, 9993964eb88SJed Brown 0, 10003964eb88SJed Brown 0, 10013964eb88SJed Brown 0, 10023964eb88SJed Brown /*139*/ 0, 1003f9426fe0SMark Adams 0, 10043964eb88SJed Brown 0 10053964eb88SJed Brown }; 1006273d9f13SBarry Smith 1007789d8953SBarry Smith PetscErrorCode MatShellSetContext_Shell(Mat mat,void *ctx) 1008789d8953SBarry Smith { 1009789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1010789d8953SBarry Smith 1011789d8953SBarry Smith PetscFunctionBegin; 1012789d8953SBarry Smith shell->ctx = ctx; 1013789d8953SBarry Smith PetscFunctionReturn(0); 1014789d8953SBarry Smith } 1015789d8953SBarry Smith 1016*db77db73SJeremy L Thompson static PetscErrorCode MatShellSetVecType_Shell(Mat mat,VecType vtype) 1017*db77db73SJeremy L Thompson { 1018*db77db73SJeremy L Thompson PetscErrorCode ierr; 1019*db77db73SJeremy L Thompson 1020*db77db73SJeremy L Thompson PetscFunctionBegin; 1021*db77db73SJeremy L Thompson ierr = PetscFree(mat->defaultvectype);CHKERRQ(ierr); 1022*db77db73SJeremy L Thompson ierr = PetscStrallocpy(vtype,(char**)&mat->defaultvectype);CHKERRQ(ierr); 1023*db77db73SJeremy L Thompson PetscFunctionReturn(0); 1024*db77db73SJeremy L Thompson } 1025*db77db73SJeremy L Thompson 1026789d8953SBarry Smith PetscErrorCode MatShellSetManageScalingShifts_Shell(Mat A) 1027789d8953SBarry Smith { 1028789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 1029789d8953SBarry Smith 1030789d8953SBarry Smith PetscFunctionBegin; 1031789d8953SBarry Smith shell->managescalingshifts = PETSC_FALSE; 1032789d8953SBarry Smith A->ops->diagonalset = NULL; 1033789d8953SBarry Smith A->ops->diagonalscale = NULL; 1034789d8953SBarry Smith A->ops->scale = NULL; 1035789d8953SBarry Smith A->ops->shift = NULL; 1036789d8953SBarry Smith A->ops->axpy = NULL; 1037789d8953SBarry Smith PetscFunctionReturn(0); 1038789d8953SBarry Smith } 1039789d8953SBarry Smith 1040789d8953SBarry Smith PetscErrorCode MatShellSetOperation_Shell(Mat mat,MatOperation op,void (*f)(void)) 1041789d8953SBarry Smith { 1042feb237baSPierre Jolivet Mat_Shell *shell = (Mat_Shell*)mat->data; 1043789d8953SBarry Smith 1044789d8953SBarry Smith PetscFunctionBegin; 1045789d8953SBarry Smith switch (op) { 1046789d8953SBarry Smith case MATOP_DESTROY: 1047789d8953SBarry Smith shell->ops->destroy = (PetscErrorCode (*)(Mat))f; 1048789d8953SBarry Smith break; 1049789d8953SBarry Smith case MATOP_VIEW: 1050789d8953SBarry Smith if (!mat->ops->viewnative) { 1051789d8953SBarry Smith mat->ops->viewnative = mat->ops->view; 1052789d8953SBarry Smith } 1053789d8953SBarry Smith mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 1054789d8953SBarry Smith break; 1055789d8953SBarry Smith case MATOP_COPY: 1056789d8953SBarry Smith shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f; 1057789d8953SBarry Smith break; 1058789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1059789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1060789d8953SBarry Smith case MATOP_SHIFT: 1061789d8953SBarry Smith case MATOP_SCALE: 1062789d8953SBarry Smith case MATOP_AXPY: 1063789d8953SBarry Smith case MATOP_ZERO_ROWS: 1064789d8953SBarry Smith case MATOP_ZERO_ROWS_COLUMNS: 1065789d8953SBarry Smith if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()"); 1066789d8953SBarry Smith (((void(**)(void))mat->ops)[op]) = f; 1067789d8953SBarry Smith break; 1068789d8953SBarry Smith case MATOP_GET_DIAGONAL: 1069789d8953SBarry Smith if (shell->managescalingshifts) { 1070789d8953SBarry Smith shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1071789d8953SBarry Smith mat->ops->getdiagonal = MatGetDiagonal_Shell; 1072789d8953SBarry Smith } else { 1073789d8953SBarry Smith shell->ops->getdiagonal = NULL; 1074789d8953SBarry Smith mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1075789d8953SBarry Smith } 1076789d8953SBarry Smith break; 1077789d8953SBarry Smith case MATOP_MULT: 1078789d8953SBarry Smith if (shell->managescalingshifts) { 1079789d8953SBarry Smith shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1080789d8953SBarry Smith mat->ops->mult = MatMult_Shell; 1081789d8953SBarry Smith } else { 1082789d8953SBarry Smith shell->ops->mult = NULL; 1083789d8953SBarry Smith mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1084789d8953SBarry Smith } 1085789d8953SBarry Smith break; 1086789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 1087789d8953SBarry Smith if (shell->managescalingshifts) { 1088789d8953SBarry Smith shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1089789d8953SBarry Smith mat->ops->multtranspose = MatMultTranspose_Shell; 1090789d8953SBarry Smith } else { 1091789d8953SBarry Smith shell->ops->multtranspose = NULL; 1092789d8953SBarry Smith mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1093789d8953SBarry Smith } 1094789d8953SBarry Smith break; 1095789d8953SBarry Smith default: 1096789d8953SBarry Smith (((void(**)(void))mat->ops)[op]) = f; 1097789d8953SBarry Smith break; 1098789d8953SBarry Smith } 1099789d8953SBarry Smith PetscFunctionReturn(0); 1100789d8953SBarry Smith } 1101789d8953SBarry Smith 1102789d8953SBarry Smith PetscErrorCode MatShellGetOperation_Shell(Mat mat,MatOperation op,void(**f)(void)) 1103789d8953SBarry Smith { 1104789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1105789d8953SBarry Smith 1106789d8953SBarry Smith PetscFunctionBegin; 1107789d8953SBarry Smith switch (op) { 1108789d8953SBarry Smith case MATOP_DESTROY: 1109789d8953SBarry Smith *f = (void (*)(void))shell->ops->destroy; 1110789d8953SBarry Smith break; 1111789d8953SBarry Smith case MATOP_VIEW: 1112789d8953SBarry Smith *f = (void (*)(void))mat->ops->view; 1113789d8953SBarry Smith break; 1114789d8953SBarry Smith case MATOP_COPY: 1115789d8953SBarry Smith *f = (void (*)(void))shell->ops->copy; 1116789d8953SBarry Smith break; 1117789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1118789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1119789d8953SBarry Smith case MATOP_SHIFT: 1120789d8953SBarry Smith case MATOP_SCALE: 1121789d8953SBarry Smith case MATOP_AXPY: 1122789d8953SBarry Smith case MATOP_ZERO_ROWS: 1123789d8953SBarry Smith case MATOP_ZERO_ROWS_COLUMNS: 1124789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1125789d8953SBarry Smith break; 1126789d8953SBarry Smith case MATOP_GET_DIAGONAL: 1127789d8953SBarry Smith if (shell->ops->getdiagonal) 1128789d8953SBarry Smith *f = (void (*)(void))shell->ops->getdiagonal; 1129789d8953SBarry Smith else 1130789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1131789d8953SBarry Smith break; 1132789d8953SBarry Smith case MATOP_MULT: 1133789d8953SBarry Smith if (shell->ops->mult) 1134789d8953SBarry Smith *f = (void (*)(void))shell->ops->mult; 1135789d8953SBarry Smith else 1136789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1137789d8953SBarry Smith break; 1138789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 1139789d8953SBarry Smith if (shell->ops->multtranspose) 1140789d8953SBarry Smith *f = (void (*)(void))shell->ops->multtranspose; 1141789d8953SBarry Smith else 1142789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1143789d8953SBarry Smith break; 1144789d8953SBarry Smith default: 1145789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1146789d8953SBarry Smith } 1147789d8953SBarry Smith PetscFunctionReturn(0); 1148789d8953SBarry Smith } 1149789d8953SBarry Smith 11500bad9183SKris Buschelman /*MC 1151fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 11520bad9183SKris Buschelman 11530bad9183SKris Buschelman Level: advanced 11540bad9183SKris Buschelman 11550c0fd78eSBarry Smith .seealso: MatCreateShell() 11560bad9183SKris Buschelman M*/ 11570bad9183SKris Buschelman 11588cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 1159273d9f13SBarry Smith { 1160273d9f13SBarry Smith Mat_Shell *b; 1161dfbe8321SBarry Smith PetscErrorCode ierr; 1162273d9f13SBarry Smith 1163273d9f13SBarry Smith PetscFunctionBegin; 1164273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1165273d9f13SBarry Smith 1166b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1167273d9f13SBarry Smith A->data = (void*)b; 1168273d9f13SBarry Smith 116926283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 117026283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1171273d9f13SBarry Smith 1172273d9f13SBarry Smith b->ctx = 0; 1173ef66eb69SBarry Smith b->vshift = 0.0; 1174ef66eb69SBarry Smith b->vscale = 1.0; 11750c0fd78eSBarry Smith b->managescalingshifts = PETSC_TRUE; 1176273d9f13SBarry Smith A->assembled = PETSC_TRUE; 1177210f0121SBarry Smith A->preallocated = PETSC_FALSE; 11782205254eSKarl Rupp 1179789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetContext_C",MatShellGetContext_Shell);CHKERRQ(ierr); 1180789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetContext_C",MatShellSetContext_Shell);CHKERRQ(ierr); 1181*db77db73SJeremy L Thompson ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetVecType_C",MatShellSetVecType_Shell);CHKERRQ(ierr); 1182789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetManageScalingShifts_C",MatShellSetManageScalingShifts_Shell);CHKERRQ(ierr); 1183789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetOperation_C",MatShellSetOperation_Shell);CHKERRQ(ierr); 1184789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetOperation_C",MatShellGetOperation_Shell);CHKERRQ(ierr); 118517667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 1186273d9f13SBarry Smith PetscFunctionReturn(0); 1187273d9f13SBarry Smith } 1188e51e0e81SBarry Smith 11894b828684SBarry Smith /*@C 1190052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 1191ff756334SLois Curfman McInnes private data storage format. 1192e51e0e81SBarry Smith 1193d083f849SBarry Smith Collective 1194c7fcc2eaSBarry Smith 1195e51e0e81SBarry Smith Input Parameters: 1196c7fcc2eaSBarry Smith + comm - MPI communicator 1197c7fcc2eaSBarry Smith . m - number of local rows (must be given) 1198c7fcc2eaSBarry Smith . n - number of local columns (must be given) 1199c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 1200c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 1201c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 1202e51e0e81SBarry Smith 1203ff756334SLois Curfman McInnes Output Parameter: 120444cd7ae7SLois Curfman McInnes . A - the matrix 1205e51e0e81SBarry Smith 1206ff2fd236SBarry Smith Level: advanced 1207ff2fd236SBarry Smith 1208f39d1f56SLois Curfman McInnes Usage: 12095bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1210f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 1211c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1212f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 1213f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 1214f39d1f56SLois Curfman McInnes 1215ff756334SLois Curfman McInnes Notes: 1216ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 1217ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 1218ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 1219e51e0e81SBarry Smith 122095452b02SPatrick Sanan Fortran Notes: 122195452b02SPatrick Sanan To use this from Fortran with a ctx you must write an interface definition for this 1222daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 1223daf670e6SBarry Smith in as the ctx argument. 1224f38a8d56SBarry Smith 1225f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 1226f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 1227645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 1228645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 1229645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 1230645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 1231645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 1232645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 1233645985a0SLois Curfman McInnes For example, 1234f39d1f56SLois Curfman McInnes 1235f39d1f56SLois Curfman McInnes $ 1236f39d1f56SLois Curfman McInnes $ Vec x, y 12375bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1238645985a0SLois Curfman McInnes $ Mat A 1239f39d1f56SLois Curfman McInnes $ 1240c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 1241c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 1242f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 1243c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 1244c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 1245c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1246645985a0SLois Curfman McInnes $ MatMult(A,x,y); 124745960306SStefano Zampini $ MatDestroy(&A); 124845960306SStefano Zampini $ VecDestroy(&y); 124945960306SStefano Zampini $ VecDestroy(&x); 1250645985a0SLois Curfman McInnes $ 1251e51e0e81SBarry Smith 12529b53d762SBarry Smith 125345960306SStefano Zampini MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), MatScale(), MatZeroRows() and MatZeroRowsColumns() internally, so these 12549b53d762SBarry Smith operations cannot be overwritten unless MatShellSetManageScalingShifts() is called. 12559b53d762SBarry Smith 12569b53d762SBarry Smith 12570c0fd78eSBarry Smith For rectangular matrices do all the scalings and shifts make sense? 12580c0fd78eSBarry Smith 125995452b02SPatrick Sanan Developers Notes: 126095452b02SPatrick Sanan Regarding shifting and scaling. The general form is 12610c0fd78eSBarry Smith 12620c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 12630c0fd78eSBarry Smith 12640c0fd78eSBarry Smith The order you apply the operations is important. For example if you have a dshift then 12650c0fd78eSBarry Smith apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift 12660c0fd78eSBarry Smith you get s*vscale*A + diag(shift) 12670c0fd78eSBarry Smith 12680c0fd78eSBarry Smith A is the user provided function. 12690c0fd78eSBarry Smith 1270ad2f5c3fSBarry Smith KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for 1271ad2f5c3fSBarry Smith for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger 1272ad2f5c3fSBarry Smith an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat); 1273ad2f5c3fSBarry Smith each time the MATSHELL matrix has changed. 1274ad2f5c3fSBarry Smith 1275ad2f5c3fSBarry Smith Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided 1276ad2f5c3fSBarry Smith with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale(). 1277ad2f5c3fSBarry Smith 12780c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts() 1279e51e0e81SBarry Smith @*/ 12807087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 1281e51e0e81SBarry Smith { 1282dfbe8321SBarry Smith PetscErrorCode ierr; 1283ed3cc1f0SBarry Smith 12843a40ed3dSBarry Smith PetscFunctionBegin; 1285f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1286f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1287273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 1288273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 12890fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 1290273d9f13SBarry Smith PetscFunctionReturn(0); 1291c7fcc2eaSBarry Smith } 1292c7fcc2eaSBarry Smith 1293789d8953SBarry Smith 1294c6866cfdSSatish Balay /*@ 1295273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 1296c7fcc2eaSBarry Smith 12973f9fe445SBarry Smith Logically Collective on Mat 1298c7fcc2eaSBarry Smith 1299273d9f13SBarry Smith Input Parameters: 1300273d9f13SBarry Smith + mat - the shell matrix 1301273d9f13SBarry Smith - ctx - the context 1302273d9f13SBarry Smith 1303273d9f13SBarry Smith Level: advanced 1304273d9f13SBarry Smith 130595452b02SPatrick Sanan Fortran Notes: 130695452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 1307daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 1308273d9f13SBarry Smith 1309273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 13100bc0a719SBarry Smith @*/ 13117087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 1312273d9f13SBarry Smith { 1313dfbe8321SBarry Smith PetscErrorCode ierr; 1314273d9f13SBarry Smith 1315273d9f13SBarry Smith PetscFunctionBegin; 13160700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1317789d8953SBarry Smith ierr = PetscUseMethod(mat,"MatShellSetContext_C",(Mat,void*),(mat,ctx));CHKERRQ(ierr); 13183a40ed3dSBarry Smith PetscFunctionReturn(0); 1319e51e0e81SBarry Smith } 1320e51e0e81SBarry Smith 1321*db77db73SJeremy L Thompson /*@C 1322*db77db73SJeremy L Thompson MatShellSetVecType - Sets the type of Vec returned by MatCreateVecs() 1323*db77db73SJeremy L Thompson 1324*db77db73SJeremy L Thompson Logically collective 1325*db77db73SJeremy L Thompson 1326*db77db73SJeremy L Thompson Input Parameters: 1327*db77db73SJeremy L Thompson + mat - the shell matrix 1328*db77db73SJeremy L Thompson - vtype - type to use for creating vectors 1329*db77db73SJeremy L Thompson 1330*db77db73SJeremy L Thompson Notes: 1331*db77db73SJeremy L Thompson 1332*db77db73SJeremy L Thompson Level: advanced 1333*db77db73SJeremy L Thompson 1334*db77db73SJeremy L Thompson .seealso: MatCreateVecs() 1335*db77db73SJeremy L Thompson @*/ 1336*db77db73SJeremy L Thompson PetscErrorCode MatShellSetVecType(Mat mat,VecType vtype) 1337*db77db73SJeremy L Thompson { 1338*db77db73SJeremy L Thompson PetscErrorCode ierr; 1339*db77db73SJeremy L Thompson 1340*db77db73SJeremy L Thompson PetscFunctionBegin; 1341*db77db73SJeremy L Thompson ierr = PetscTryMethod(mat,"MatShellSetVecType_C",(Mat,VecType),(mat,vtype));CHKERRQ(ierr); 1342*db77db73SJeremy L Thompson PetscFunctionReturn(0); 1343*db77db73SJeremy L Thompson } 1344*db77db73SJeremy L Thompson 13450c0fd78eSBarry Smith /*@ 13460c0fd78eSBarry Smith MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately 13470c0fd78eSBarry Smith after MatCreateShell() 13480c0fd78eSBarry Smith 13490c0fd78eSBarry Smith Logically Collective on Mat 13500c0fd78eSBarry Smith 13510c0fd78eSBarry Smith Input Parameter: 13520c0fd78eSBarry Smith . mat - the shell matrix 13530c0fd78eSBarry Smith 13540c0fd78eSBarry Smith Level: advanced 13550c0fd78eSBarry Smith 13560c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation() 13570c0fd78eSBarry Smith @*/ 13580c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A) 13590c0fd78eSBarry Smith { 13600c0fd78eSBarry Smith PetscErrorCode ierr; 13610c0fd78eSBarry Smith 13620c0fd78eSBarry Smith PetscFunctionBegin; 13630c0fd78eSBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1364789d8953SBarry Smith ierr = PetscUseMethod(A,"MatShellSetManageScalingShifts_C",(Mat),(A));CHKERRQ(ierr); 13650c0fd78eSBarry Smith PetscFunctionReturn(0); 13660c0fd78eSBarry Smith } 13670c0fd78eSBarry Smith 1368c16cb8f2SBarry Smith /*@C 1369f3b1f45cSBarry Smith MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function. 1370f3b1f45cSBarry Smith 1371f3b1f45cSBarry Smith Logically Collective on Mat 1372f3b1f45cSBarry Smith 1373f3b1f45cSBarry Smith Input Parameters: 1374f3b1f45cSBarry Smith + mat - the shell matrix 1375f3b1f45cSBarry Smith . f - the function 1376f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1377f3b1f45cSBarry Smith - ctx - an optional context for the function 1378f3b1f45cSBarry Smith 1379f3b1f45cSBarry Smith Output Parameter: 1380f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1381f3b1f45cSBarry Smith 1382f3b1f45cSBarry Smith Options Database: 1383f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1384f3b1f45cSBarry Smith 1385f3b1f45cSBarry Smith Level: advanced 1386f3b1f45cSBarry Smith 138795452b02SPatrick Sanan Fortran Notes: 138895452b02SPatrick Sanan Not supported from Fortran 1389f3b1f45cSBarry Smith 1390f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose() 1391f3b1f45cSBarry Smith @*/ 1392f3b1f45cSBarry Smith PetscErrorCode MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1393f3b1f45cSBarry Smith { 1394f3b1f45cSBarry Smith PetscErrorCode ierr; 1395f3b1f45cSBarry Smith PetscInt m,n; 1396f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1397f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 139874e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1399f3b1f45cSBarry Smith 1400f3b1f45cSBarry Smith PetscFunctionBegin; 1401f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1402f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr); 1403f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 1404f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr); 1405f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 1406f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 1407f3b1f45cSBarry Smith 14080bacdadaSStefano Zampini ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr); 14090bacdadaSStefano Zampini ierr = MatComputeOperator(mat,MATAIJ,&Dmat);CHKERRQ(ierr); 1410f3b1f45cSBarry Smith 1411f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 1412f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 1413f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 1414f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 1415f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1416f3b1f45cSBarry Smith flag = PETSC_FALSE; 1417f3b1f45cSBarry Smith if (v) { 1418fc7aafd1SBarry 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); 1419f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 1420f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 1421f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 1422f3b1f45cSBarry Smith } 1423f3b1f45cSBarry Smith } else if (v) { 1424fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr); 1425f3b1f45cSBarry Smith } 1426f3b1f45cSBarry Smith if (flg) *flg = flag; 1427f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 1428f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 1429f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 1430f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 1431f3b1f45cSBarry Smith PetscFunctionReturn(0); 1432f3b1f45cSBarry Smith } 1433f3b1f45cSBarry Smith 1434f3b1f45cSBarry Smith /*@C 1435f3b1f45cSBarry Smith MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function. 1436f3b1f45cSBarry Smith 1437f3b1f45cSBarry Smith Logically Collective on Mat 1438f3b1f45cSBarry Smith 1439f3b1f45cSBarry Smith Input Parameters: 1440f3b1f45cSBarry Smith + mat - the shell matrix 1441f3b1f45cSBarry Smith . f - the function 1442f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1443f3b1f45cSBarry Smith - ctx - an optional context for the function 1444f3b1f45cSBarry Smith 1445f3b1f45cSBarry Smith Output Parameter: 1446f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1447f3b1f45cSBarry Smith 1448f3b1f45cSBarry Smith Options Database: 1449f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1450f3b1f45cSBarry Smith 1451f3b1f45cSBarry Smith Level: advanced 1452f3b1f45cSBarry Smith 145395452b02SPatrick Sanan Fortran Notes: 145495452b02SPatrick Sanan Not supported from Fortran 1455f3b1f45cSBarry Smith 1456f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult() 1457f3b1f45cSBarry Smith @*/ 1458f3b1f45cSBarry Smith PetscErrorCode MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1459f3b1f45cSBarry Smith { 1460f3b1f45cSBarry Smith PetscErrorCode ierr; 1461f3b1f45cSBarry Smith Vec x,y,z; 1462f3b1f45cSBarry Smith PetscInt m,n,M,N; 1463f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1464f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 146574e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1466f3b1f45cSBarry Smith 1467f3b1f45cSBarry Smith PetscFunctionBegin; 1468f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1469f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr); 1470f3b1f45cSBarry Smith ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr); 1471f3b1f45cSBarry Smith ierr = VecDuplicate(y,&z);CHKERRQ(ierr); 1472f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 1473f3b1f45cSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 1474f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr); 1475f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 1476f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 14770bacdadaSStefano Zampini ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr); 1478f3b1f45cSBarry Smith ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr); 14790bacdadaSStefano Zampini ierr = MatComputeOperatorTranspose(mat,MATAIJ,&Dmat);CHKERRQ(ierr); 1480f3b1f45cSBarry Smith 1481f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 1482f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 1483f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 1484f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 1485f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1486f3b1f45cSBarry Smith flag = PETSC_FALSE; 1487f3b1f45cSBarry Smith if (v) { 1488fc7aafd1SBarry 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); 1489f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 1490f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 1491f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 1492f3b1f45cSBarry Smith } 1493f3b1f45cSBarry Smith } else if (v) { 1494fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr); 1495f3b1f45cSBarry Smith } 1496f3b1f45cSBarry Smith if (flg) *flg = flag; 1497f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 1498f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 1499f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 1500f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 1501f3b1f45cSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1502f3b1f45cSBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 1503f3b1f45cSBarry Smith ierr = VecDestroy(&z);CHKERRQ(ierr); 1504f3b1f45cSBarry Smith PetscFunctionReturn(0); 1505f3b1f45cSBarry Smith } 1506f3b1f45cSBarry Smith 1507f3b1f45cSBarry Smith /*@C 15080c0fd78eSBarry Smith MatShellSetOperation - Allows user to set a matrix operation for a shell matrix. 1509e51e0e81SBarry Smith 15103f9fe445SBarry Smith Logically Collective on Mat 1511fee21e36SBarry Smith 1512c7fcc2eaSBarry Smith Input Parameters: 1513c7fcc2eaSBarry Smith + mat - the shell matrix 1514c7fcc2eaSBarry Smith . op - the name of the operation 1515789d8953SBarry Smith - g - the function that provides the operation. 1516c7fcc2eaSBarry Smith 151715091d37SBarry Smith Level: advanced 151815091d37SBarry Smith 1519fae171e0SBarry Smith Usage: 1520e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 1521f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 1522c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 15230b627109SLois Curfman McInnes 1524a62d957aSLois Curfman McInnes Notes: 1525e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 15261c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 1527a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 15281c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 1529a62d957aSLois Curfman McInnes 1530a4d85fd7SBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 1531deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 1532deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 1533deebb3c3SLois Curfman McInnes routines, e.g., 1534a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1535a62d957aSLois Curfman McInnes 15364aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 15374aa34b0aSBarry Smith nonzero on failure. 15384aa34b0aSBarry Smith 1539a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 1540a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 1541a62d957aSLois Curfman McInnes set by MatCreateShell(). 1542a62d957aSLois Curfman McInnes 154395452b02SPatrick Sanan Fortran Notes: 154495452b02SPatrick Sanan For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 1545501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 1546501d9185SBarry Smith 15470c0fd78eSBarry Smith Use MatSetOperation() to set an operation for any matrix type 15480c0fd78eSBarry Smith 15490c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts() 1550e51e0e81SBarry Smith @*/ 1551789d8953SBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*g)(void)) 1552e51e0e81SBarry Smith { 1553976e8c5aSLisandro Dalcin PetscErrorCode ierr; 1554273d9f13SBarry Smith 15553a40ed3dSBarry Smith PetscFunctionBegin; 15560700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1557789d8953SBarry Smith ierr = PetscUseMethod(mat,"MatShellSetOperation_C",(Mat,MatOperation,void (*)(void)),(mat,op,g));CHKERRQ(ierr); 15583a40ed3dSBarry Smith PetscFunctionReturn(0); 1559e51e0e81SBarry Smith } 1560f0479e8cSBarry Smith 1561d4bb536fSBarry Smith /*@C 1562d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 1563d4bb536fSBarry Smith 1564c7fcc2eaSBarry Smith Not Collective 1565c7fcc2eaSBarry Smith 1566d4bb536fSBarry Smith Input Parameters: 1567c7fcc2eaSBarry Smith + mat - the shell matrix 1568c7fcc2eaSBarry Smith - op - the name of the operation 1569d4bb536fSBarry Smith 1570d4bb536fSBarry Smith Output Parameter: 1571789d8953SBarry Smith . g - the function that provides the operation. 1572d4bb536fSBarry Smith 157315091d37SBarry Smith Level: advanced 157415091d37SBarry Smith 1575d4bb536fSBarry Smith Notes: 1576e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 1577d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 1578d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 1579d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 1580d4bb536fSBarry Smith 1581d4bb536fSBarry Smith All user-provided functions have the same calling 1582d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 1583d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 1584d4bb536fSBarry Smith routines, e.g., 1585d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1586d4bb536fSBarry Smith 1587d4bb536fSBarry Smith Within each user-defined routine, the user should call 1588d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 1589d4bb536fSBarry Smith set by MatCreateShell(). 1590d4bb536fSBarry Smith 1591ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 1592d4bb536fSBarry Smith @*/ 1593789d8953SBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**g)(void)) 1594d4bb536fSBarry Smith { 1595976e8c5aSLisandro Dalcin PetscErrorCode ierr; 1596273d9f13SBarry Smith 15973a40ed3dSBarry Smith PetscFunctionBegin; 15980700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1599789d8953SBarry Smith ierr = PetscUseMethod(mat,"MatShellGetOperation_C",(Mat,MatOperation,void (**)(void)),(mat,op,g));CHKERRQ(ierr); 16003a40ed3dSBarry Smith PetscFunctionReturn(0); 1601d4bb536fSBarry Smith } 1602