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); 476789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetManageScalingShifts_C",NULL);CHKERRQ(ierr); 477789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetOperation_C",NULL);CHKERRQ(ierr); 478789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellGetOperation_C",NULL);CHKERRQ(ierr); 4793a40ed3dSBarry Smith PetscFunctionReturn(0); 480e51e0e81SBarry Smith } 481e51e0e81SBarry Smith 4827fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str) 4837fabda3fSAlex Fikl { 48428f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data; 4857fabda3fSAlex Fikl PetscErrorCode ierr; 486cb8c8a77SBarry Smith PetscBool matflg; 4877fabda3fSAlex Fikl 4887fabda3fSAlex Fikl PetscFunctionBegin; 48928f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr); 490cb8c8a77SBarry Smith if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL"); 4917fabda3fSAlex Fikl 492cb8c8a77SBarry Smith ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 493cb8c8a77SBarry Smith ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr); 4947fabda3fSAlex Fikl 495cb8c8a77SBarry Smith if (shellA->ops->copy) { 49628f357e3SAlex Fikl ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr); 497cb8c8a77SBarry Smith } 4987fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 4997fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 5000c0fd78eSBarry Smith if (shellA->dshift) { 5010c0fd78eSBarry Smith if (!shellB->dshift) { 5020c0fd78eSBarry Smith ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr); 5037fabda3fSAlex Fikl } 5040c0fd78eSBarry Smith ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr); 5057fabda3fSAlex Fikl } else { 5069f137db4SBarry Smith ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr); 5077fabda3fSAlex Fikl } 5080c0fd78eSBarry Smith if (shellA->left) { 5090c0fd78eSBarry Smith if (!shellB->left) { 5100c0fd78eSBarry Smith ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr); 5117fabda3fSAlex Fikl } 5120c0fd78eSBarry Smith ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr); 5137fabda3fSAlex Fikl } else { 5149f137db4SBarry Smith ierr = VecDestroy(&shellB->left);CHKERRQ(ierr); 5157fabda3fSAlex Fikl } 5160c0fd78eSBarry Smith if (shellA->right) { 5170c0fd78eSBarry Smith if (!shellB->right) { 5180c0fd78eSBarry Smith ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr); 5197fabda3fSAlex Fikl } 5200c0fd78eSBarry Smith ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr); 5217fabda3fSAlex Fikl } else { 5229f137db4SBarry Smith ierr = VecDestroy(&shellB->right);CHKERRQ(ierr); 5237fabda3fSAlex Fikl } 52440e381c3SBarry Smith ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr); 52540e381c3SBarry Smith if (shellA->axpy) { 52640e381c3SBarry Smith ierr = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr); 52740e381c3SBarry Smith shellB->axpy = shellA->axpy; 52840e381c3SBarry Smith shellB->axpy_vscale = shellA->axpy_vscale; 52940e381c3SBarry Smith } 53045960306SStefano Zampini if (shellA->zrows) { 53145960306SStefano Zampini ierr = ISDuplicate(shellA->zrows,&shellB->zrows);CHKERRQ(ierr); 53245960306SStefano Zampini if (shellA->zcols) { 53345960306SStefano Zampini ierr = ISDuplicate(shellA->zcols,&shellB->zcols);CHKERRQ(ierr); 53445960306SStefano Zampini } 53545960306SStefano Zampini ierr = VecDuplicate(shellA->zvals,&shellB->zvals);CHKERRQ(ierr); 53645960306SStefano Zampini ierr = VecCopy(shellA->zvals,shellB->zvals);CHKERRQ(ierr); 53745960306SStefano Zampini ierr = VecDuplicate(shellA->zvals_w,&shellB->zvals_w);CHKERRQ(ierr); 53845960306SStefano Zampini ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_r);CHKERRQ(ierr); 53945960306SStefano Zampini ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_c);CHKERRQ(ierr); 54045960306SStefano Zampini shellB->zvals_sct_r = shellA->zvals_sct_r; 54145960306SStefano Zampini shellB->zvals_sct_c = shellA->zvals_sct_c; 54245960306SStefano Zampini } 5437fabda3fSAlex Fikl PetscFunctionReturn(0); 5447fabda3fSAlex Fikl } 5457fabda3fSAlex Fikl 546cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M) 547cb8c8a77SBarry Smith { 548cb8c8a77SBarry Smith PetscErrorCode ierr; 549cb8c8a77SBarry Smith void *ctx; 550cb8c8a77SBarry Smith 551cb8c8a77SBarry Smith PetscFunctionBegin; 552cb8c8a77SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 553cb8c8a77SBarry Smith ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr); 554a4b1380bSStefano Zampini if (op != MAT_DO_NOT_COPY_VALUES) { 555cb8c8a77SBarry Smith ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 556a4b1380bSStefano Zampini } 557cb8c8a77SBarry Smith PetscFunctionReturn(0); 558cb8c8a77SBarry Smith } 559cb8c8a77SBarry Smith 560dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 561ef66eb69SBarry Smith { 562ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 563dfbe8321SBarry Smith PetscErrorCode ierr; 56425578ef6SJed Brown Vec xx; 565e3f487b0SBarry Smith PetscObjectState instate,outstate; 566ef66eb69SBarry Smith 567ef66eb69SBarry Smith PetscFunctionBegin; 568976e8c5aSLisandro Dalcin if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL"); 56945960306SStefano Zampini ierr = MatShellPreZeroRight(A,x,&xx);CHKERRQ(ierr); 57045960306SStefano Zampini ierr = MatShellPreScaleRight(A,xx,&xx);CHKERRQ(ierr); 57145960306SStefano Zampini ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 57228f357e3SAlex Fikl ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr); 573e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 574e3f487b0SBarry Smith if (instate == outstate) { 575e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 576e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 577e3f487b0SBarry Smith } 5788fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 5798fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 58045960306SStefano Zampini ierr = MatShellPostZeroLeft(A,y);CHKERRQ(ierr); 5819f137db4SBarry Smith 5829f137db4SBarry Smith if (shell->axpy) { 5839f137db4SBarry Smith if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);} 5849f137db4SBarry Smith ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr); 5859f137db4SBarry Smith ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 5869f137db4SBarry Smith } 587ef66eb69SBarry Smith PetscFunctionReturn(0); 588ef66eb69SBarry Smith } 589ef66eb69SBarry Smith 5905edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 5915edf6516SJed Brown { 5925edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 5935edf6516SJed Brown PetscErrorCode ierr; 5945edf6516SJed Brown 5955edf6516SJed Brown PetscFunctionBegin; 5965edf6516SJed Brown if (y == z) { 5975edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 5985edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 599b8430d49SBarry Smith ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr); 6005edf6516SJed Brown } else { 6015edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 6025edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 6035edf6516SJed Brown } 6045edf6516SJed Brown PetscFunctionReturn(0); 6055edf6516SJed Brown } 6065edf6516SJed Brown 60718be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 60818be62a5SSatish Balay { 60918be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 61018be62a5SSatish Balay PetscErrorCode ierr; 61125578ef6SJed Brown Vec xx; 612e3f487b0SBarry Smith PetscObjectState instate,outstate; 61318be62a5SSatish Balay 61418be62a5SSatish Balay PetscFunctionBegin; 61545960306SStefano Zampini ierr = MatShellPreZeroLeft(A,x,&xx);CHKERRQ(ierr); 61645960306SStefano Zampini ierr = MatShellPreScaleLeft(A,xx,&xx);CHKERRQ(ierr); 617e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 6180c0fd78eSBarry Smith if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL"); 61928f357e3SAlex Fikl ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr); 620e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 621e3f487b0SBarry Smith if (instate == outstate) { 622e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 623e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 624e3f487b0SBarry Smith } 6258fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 6268fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 62745960306SStefano Zampini ierr = MatShellPostZeroRight(A,y);CHKERRQ(ierr); 628350ec83bSStefano Zampini 629350ec83bSStefano Zampini if (shell->axpy) { 630350ec83bSStefano Zampini if (!shell->right_work) {ierr = MatCreateVecs(A,NULL,&shell->right_work);CHKERRQ(ierr);} 631350ec83bSStefano Zampini ierr = MatMultTranspose(shell->axpy,x,shell->right_work);CHKERRQ(ierr); 632350ec83bSStefano Zampini ierr = VecAXPY(y,shell->axpy_vscale,shell->right_work);CHKERRQ(ierr); 633350ec83bSStefano Zampini } 63418be62a5SSatish Balay PetscFunctionReturn(0); 63518be62a5SSatish Balay } 63618be62a5SSatish Balay 6375edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 6385edf6516SJed Brown { 6395edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 6405edf6516SJed Brown PetscErrorCode ierr; 6415edf6516SJed Brown 6425edf6516SJed Brown PetscFunctionBegin; 6435edf6516SJed Brown if (y == z) { 6445edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 6455edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 646dac989eaS“Marek ierr = VecAXPY(z,1.0,shell->left_add_work);CHKERRQ(ierr); 6475edf6516SJed Brown } else { 6485edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 6495edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 6505edf6516SJed Brown } 6515edf6516SJed Brown PetscFunctionReturn(0); 6525edf6516SJed Brown } 6535edf6516SJed Brown 6540c0fd78eSBarry Smith /* 6550c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 6560c0fd78eSBarry Smith */ 65718be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 65818be62a5SSatish Balay { 65918be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 66018be62a5SSatish Balay PetscErrorCode ierr; 66118be62a5SSatish Balay 66218be62a5SSatish Balay PetscFunctionBegin; 6630c0fd78eSBarry Smith if (shell->ops->getdiagonal) { 66428f357e3SAlex Fikl ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr); 665305211d5SBarry 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,...)"); 66618be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 6678fe8eb27SJed Brown if (shell->dshift) { 6680c0fd78eSBarry Smith ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr); 66918be62a5SSatish Balay } 6700c0fd78eSBarry Smith ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 6718fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 6728fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 67345960306SStefano Zampini if (shell->zrows) { 67445960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 67545960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 67645960306SStefano Zampini } 67781c02519SBarry Smith if (shell->axpy) { 67881c02519SBarry Smith if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);} 67981c02519SBarry Smith ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr); 68081c02519SBarry Smith ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 68181c02519SBarry Smith } 68218be62a5SSatish Balay PetscFunctionReturn(0); 68318be62a5SSatish Balay } 68418be62a5SSatish Balay 685f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 686ef66eb69SBarry Smith { 687ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 6888fe8eb27SJed Brown PetscErrorCode ierr; 689789d8953SBarry Smith PetscBool flg; 690b24ad042SBarry Smith 691ef66eb69SBarry Smith PetscFunctionBegin; 692789d8953SBarry Smith ierr = MatHasCongruentLayouts(Y,&flg);CHKERRQ(ierr); 693789d8953SBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot shift shell matrix if it is not congruent"); 6940c0fd78eSBarry Smith if (shell->left || shell->right) { 6958fe8eb27SJed Brown if (!shell->dshift) { 6960c0fd78eSBarry Smith ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr); 6970c0fd78eSBarry Smith ierr = VecSet(shell->dshift,a);CHKERRQ(ierr); 6980c0fd78eSBarry Smith } else { 6990c0fd78eSBarry Smith if (shell->left) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 7000c0fd78eSBarry Smith if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 7019f137db4SBarry Smith ierr = VecShift(shell->dshift,a);CHKERRQ(ierr); 7020c0fd78eSBarry Smith } 7038fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 7048fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 7058fe8eb27SJed Brown } else shell->vshift += a; 70645960306SStefano Zampini if (shell->zrows) { 70745960306SStefano Zampini ierr = VecShift(shell->zvals,a);CHKERRQ(ierr); 70845960306SStefano Zampini } 709ef66eb69SBarry Smith PetscFunctionReturn(0); 710ef66eb69SBarry Smith } 711ef66eb69SBarry Smith 712b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s) 7136464bf51SAlex Fikl { 7146464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 7156464bf51SAlex Fikl PetscErrorCode ierr; 7166464bf51SAlex Fikl 7176464bf51SAlex Fikl PetscFunctionBegin; 7180c0fd78eSBarry Smith if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);} 7190c0fd78eSBarry Smith if (shell->left || shell->right) { 72092fabd1fSBarry Smith if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);} 7219f137db4SBarry Smith if (shell->left && shell->right) { 7229f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 7239f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr); 7249f137db4SBarry Smith } else if (shell->left) { 7259f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 7269f137db4SBarry Smith } else { 7279f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr); 7289f137db4SBarry Smith } 729b253ae0bS“Marek ierr = VecAXPY(shell->dshift,s,shell->right_work);CHKERRQ(ierr); 7300c0fd78eSBarry Smith } else { 731b253ae0bS“Marek ierr = VecAXPY(shell->dshift,s,D);CHKERRQ(ierr); 732b253ae0bS“Marek } 733b253ae0bS“Marek PetscFunctionReturn(0); 734b253ae0bS“Marek } 735b253ae0bS“Marek 736b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 737b253ae0bS“Marek { 73845960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 739b253ae0bS“Marek Vec d; 740b253ae0bS“Marek PetscErrorCode ierr; 741789d8953SBarry Smith PetscBool flg; 742b253ae0bS“Marek 743b253ae0bS“Marek PetscFunctionBegin; 744789d8953SBarry Smith ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr); 745789d8953SBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot diagonal set or shift shell matrix if it is not congruent"); 746b253ae0bS“Marek if (ins == INSERT_VALUES) { 747b253ae0bS“Marek if (!A->ops->getdiagonal) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first"); 748b253ae0bS“Marek ierr = VecDuplicate(D,&d);CHKERRQ(ierr); 749b253ae0bS“Marek ierr = MatGetDiagonal(A,d);CHKERRQ(ierr); 750b253ae0bS“Marek ierr = MatDiagonalSet_Shell_Private(A,d,-1.);CHKERRQ(ierr); 751b253ae0bS“Marek ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr); 752b253ae0bS“Marek ierr = VecDestroy(&d);CHKERRQ(ierr); 75345960306SStefano Zampini if (shell->zrows) { 75445960306SStefano Zampini ierr = VecCopy(D,shell->zvals);CHKERRQ(ierr); 75545960306SStefano Zampini } 756b253ae0bS“Marek } else { 757b253ae0bS“Marek ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr); 75845960306SStefano Zampini if (shell->zrows) { 75945960306SStefano Zampini ierr = VecAXPY(shell->zvals,1.0,D);CHKERRQ(ierr); 76045960306SStefano Zampini } 7616464bf51SAlex Fikl } 7626464bf51SAlex Fikl PetscFunctionReturn(0); 7636464bf51SAlex Fikl } 7646464bf51SAlex Fikl 765f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 766ef66eb69SBarry Smith { 767ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 7688fe8eb27SJed Brown PetscErrorCode ierr; 769b24ad042SBarry Smith 770ef66eb69SBarry Smith PetscFunctionBegin; 771f4df32b1SMatthew Knepley shell->vscale *= a; 7720c0fd78eSBarry Smith shell->vshift *= a; 7732205254eSKarl Rupp if (shell->dshift) { 7742205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 7750c0fd78eSBarry Smith } 77681c02519SBarry Smith shell->axpy_vscale *= a; 77745960306SStefano Zampini if (shell->zrows) { 77845960306SStefano Zampini ierr = VecScale(shell->zvals,a);CHKERRQ(ierr); 77945960306SStefano Zampini } 7808fe8eb27SJed Brown PetscFunctionReturn(0); 78118be62a5SSatish Balay } 7828fe8eb27SJed Brown 7838fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 7848fe8eb27SJed Brown { 7858fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 7868fe8eb27SJed Brown PetscErrorCode ierr; 7878fe8eb27SJed Brown 7888fe8eb27SJed Brown PetscFunctionBegin; 78981c02519SBarry Smith if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation"); 7908fe8eb27SJed Brown if (left) { 7910c0fd78eSBarry Smith if (!shell->left) { 7920c0fd78eSBarry Smith ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr); 7938fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 7940c0fd78eSBarry Smith } else { 7950c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 79618be62a5SSatish Balay } 79745960306SStefano Zampini if (shell->zrows) { 79845960306SStefano Zampini ierr = VecPointwiseMult(shell->zvals,shell->zvals,left);CHKERRQ(ierr); 79945960306SStefano Zampini } 800ef66eb69SBarry Smith } 8018fe8eb27SJed Brown if (right) { 8020c0fd78eSBarry Smith if (!shell->right) { 8030c0fd78eSBarry Smith ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr); 8048fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 8050c0fd78eSBarry Smith } else { 8060c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 8078fe8eb27SJed Brown } 80845960306SStefano Zampini if (shell->zrows) { 80945960306SStefano Zampini if (!shell->left_work) { 81045960306SStefano Zampini ierr = MatCreateVecs(Y,NULL,&shell->left_work);CHKERRQ(ierr); 81145960306SStefano Zampini } 81245960306SStefano Zampini ierr = VecSet(shell->zvals_w,1.0);CHKERRQ(ierr); 81345960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 81445960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 81545960306SStefano Zampini ierr = VecPointwiseMult(shell->zvals,shell->zvals,shell->zvals_w);CHKERRQ(ierr); 81645960306SStefano Zampini } 8178fe8eb27SJed Brown } 818ef66eb69SBarry Smith PetscFunctionReturn(0); 819ef66eb69SBarry Smith } 820ef66eb69SBarry Smith 821dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 822ef66eb69SBarry Smith { 823ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 8240c0fd78eSBarry Smith PetscErrorCode ierr; 825ef66eb69SBarry Smith 826ef66eb69SBarry Smith PetscFunctionBegin; 8278fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 828ef66eb69SBarry Smith shell->vshift = 0.0; 829ef66eb69SBarry Smith shell->vscale = 1.0; 8300c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 8310c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 8320c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 83381c02519SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 83445960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr); 83545960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr); 83645960306SStefano Zampini ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr); 83745960306SStefano Zampini ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr); 838ef66eb69SBarry Smith } 839ef66eb69SBarry Smith PetscFunctionReturn(0); 840ef66eb69SBarry Smith } 841ef66eb69SBarry Smith 8423b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 8433b49f96aSBarry Smith { 8443b49f96aSBarry Smith PetscFunctionBegin; 8453b49f96aSBarry Smith *missing = PETSC_FALSE; 8463b49f96aSBarry Smith PetscFunctionReturn(0); 8473b49f96aSBarry Smith } 8483b49f96aSBarry Smith 8499f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str) 8509f137db4SBarry Smith { 8519f137db4SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 8529f137db4SBarry Smith PetscErrorCode ierr; 8539f137db4SBarry Smith 8549f137db4SBarry Smith PetscFunctionBegin; 8559f137db4SBarry Smith ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 8569f137db4SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 8579f137db4SBarry Smith shell->axpy = X; 8589f137db4SBarry Smith shell->axpy_vscale = a; 8599f137db4SBarry Smith PetscFunctionReturn(0); 8609f137db4SBarry Smith } 8619f137db4SBarry Smith 86209dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 86320563c6bSBarry Smith 0, 86420563c6bSBarry Smith 0, 8659f137db4SBarry Smith 0, 8660c0fd78eSBarry Smith /* 4*/ MatMultAdd_Shell, 8679f137db4SBarry Smith 0, 8680c0fd78eSBarry Smith MatMultTransposeAdd_Shell, 869b951964fSBarry Smith 0, 870b951964fSBarry Smith 0, 871b951964fSBarry Smith 0, 87297304618SKris Buschelman /*10*/ 0, 873b951964fSBarry Smith 0, 874b951964fSBarry Smith 0, 875b951964fSBarry Smith 0, 876b951964fSBarry Smith 0, 87797304618SKris Buschelman /*15*/ 0, 878b951964fSBarry Smith 0, 87900849d43SBarry Smith 0, 8808fe8eb27SJed Brown MatDiagonalScale_Shell, 881b951964fSBarry Smith 0, 88297304618SKris Buschelman /*20*/ 0, 883ef66eb69SBarry Smith MatAssemblyEnd_Shell, 884b951964fSBarry Smith 0, 885b951964fSBarry Smith 0, 88645960306SStefano Zampini /*24*/ MatZeroRows_Shell, 887b951964fSBarry Smith 0, 888b951964fSBarry Smith 0, 889b951964fSBarry Smith 0, 890b951964fSBarry Smith 0, 891d519adbfSMatthew Knepley /*29*/ 0, 892b951964fSBarry Smith 0, 893273d9f13SBarry Smith 0, 894b951964fSBarry Smith 0, 895b951964fSBarry Smith 0, 896cb8c8a77SBarry Smith /*34*/ MatDuplicate_Shell, 897b951964fSBarry Smith 0, 898b951964fSBarry Smith 0, 89909dc0095SBarry Smith 0, 90009dc0095SBarry Smith 0, 9019f137db4SBarry Smith /*39*/ MatAXPY_Shell, 90209dc0095SBarry Smith 0, 90309dc0095SBarry Smith 0, 90409dc0095SBarry Smith 0, 905cb8c8a77SBarry Smith MatCopy_Shell, 906d519adbfSMatthew Knepley /*44*/ 0, 907ef66eb69SBarry Smith MatScale_Shell, 908ef66eb69SBarry Smith MatShift_Shell, 9096464bf51SAlex Fikl MatDiagonalSet_Shell, 91045960306SStefano Zampini MatZeroRowsColumns_Shell, 911f73d5cc4SBarry Smith /*49*/ 0, 91209dc0095SBarry Smith 0, 91309dc0095SBarry Smith 0, 91409dc0095SBarry Smith 0, 91509dc0095SBarry Smith 0, 916d519adbfSMatthew Knepley /*54*/ 0, 91709dc0095SBarry Smith 0, 91809dc0095SBarry Smith 0, 91909dc0095SBarry Smith 0, 92009dc0095SBarry Smith 0, 921d519adbfSMatthew Knepley /*59*/ 0, 922b9b97703SBarry Smith MatDestroy_Shell, 92309dc0095SBarry Smith 0, 924251fad3fSStefano Zampini MatConvertFrom_Shell, 925273d9f13SBarry Smith 0, 926d519adbfSMatthew Knepley /*64*/ 0, 927273d9f13SBarry Smith 0, 928273d9f13SBarry Smith 0, 929273d9f13SBarry Smith 0, 930273d9f13SBarry Smith 0, 931d519adbfSMatthew Knepley /*69*/ 0, 932273d9f13SBarry Smith 0, 933c87e5d42SMatthew Knepley MatConvert_Shell, 934273d9f13SBarry Smith 0, 93597304618SKris Buschelman 0, 936d519adbfSMatthew Knepley /*74*/ 0, 93797304618SKris Buschelman 0, 93897304618SKris Buschelman 0, 93997304618SKris Buschelman 0, 94097304618SKris Buschelman 0, 941d519adbfSMatthew Knepley /*79*/ 0, 94297304618SKris Buschelman 0, 94397304618SKris Buschelman 0, 94497304618SKris Buschelman 0, 945865e5f61SKris Buschelman 0, 946d519adbfSMatthew Knepley /*84*/ 0, 947865e5f61SKris Buschelman 0, 948865e5f61SKris Buschelman 0, 949865e5f61SKris Buschelman 0, 950865e5f61SKris Buschelman 0, 951d519adbfSMatthew Knepley /*89*/ 0, 952865e5f61SKris Buschelman 0, 953865e5f61SKris Buschelman 0, 954865e5f61SKris Buschelman 0, 955865e5f61SKris Buschelman 0, 956d519adbfSMatthew Knepley /*94*/ 0, 957865e5f61SKris Buschelman 0, 958865e5f61SKris Buschelman 0, 9593964eb88SJed Brown 0, 9603964eb88SJed Brown 0, 9613964eb88SJed Brown /*99*/ 0, 9623964eb88SJed Brown 0, 9633964eb88SJed Brown 0, 9643964eb88SJed Brown 0, 9653964eb88SJed Brown 0, 9663964eb88SJed Brown /*104*/ 0, 9673964eb88SJed Brown 0, 9683964eb88SJed Brown 0, 9693964eb88SJed Brown 0, 9703964eb88SJed Brown 0, 9713964eb88SJed Brown /*109*/ 0, 9723964eb88SJed Brown 0, 9733964eb88SJed Brown 0, 9743964eb88SJed Brown 0, 9753b49f96aSBarry Smith MatMissingDiagonal_Shell, 9763964eb88SJed Brown /*114*/ 0, 9773964eb88SJed Brown 0, 9783964eb88SJed Brown 0, 9793964eb88SJed Brown 0, 9803964eb88SJed Brown 0, 9813964eb88SJed Brown /*119*/ 0, 9823964eb88SJed Brown 0, 9833964eb88SJed Brown 0, 9843964eb88SJed Brown 0, 9853964eb88SJed Brown 0, 9863964eb88SJed Brown /*124*/ 0, 9873964eb88SJed Brown 0, 9883964eb88SJed Brown 0, 9893964eb88SJed Brown 0, 9903964eb88SJed Brown 0, 9913964eb88SJed Brown /*129*/ 0, 9923964eb88SJed Brown 0, 9933964eb88SJed Brown 0, 9943964eb88SJed Brown 0, 9953964eb88SJed Brown 0, 9963964eb88SJed Brown /*134*/ 0, 9973964eb88SJed Brown 0, 9983964eb88SJed Brown 0, 9993964eb88SJed Brown 0, 10003964eb88SJed Brown 0, 10013964eb88SJed Brown /*139*/ 0, 1002f9426fe0SMark Adams 0, 10033964eb88SJed Brown 0 10043964eb88SJed Brown }; 1005273d9f13SBarry Smith 1006789d8953SBarry Smith PetscErrorCode MatShellSetContext_Shell(Mat mat,void *ctx) 1007789d8953SBarry Smith { 1008789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1009789d8953SBarry Smith 1010789d8953SBarry Smith PetscFunctionBegin; 1011789d8953SBarry Smith shell->ctx = ctx; 1012789d8953SBarry Smith PetscFunctionReturn(0); 1013789d8953SBarry Smith } 1014789d8953SBarry Smith 1015789d8953SBarry Smith PetscErrorCode MatShellSetManageScalingShifts_Shell(Mat A) 1016789d8953SBarry Smith { 1017789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 1018789d8953SBarry Smith 1019789d8953SBarry Smith PetscFunctionBegin; 1020789d8953SBarry Smith shell->managescalingshifts = PETSC_FALSE; 1021789d8953SBarry Smith A->ops->diagonalset = NULL; 1022789d8953SBarry Smith A->ops->diagonalscale = NULL; 1023789d8953SBarry Smith A->ops->scale = NULL; 1024789d8953SBarry Smith A->ops->shift = NULL; 1025789d8953SBarry Smith A->ops->axpy = NULL; 1026789d8953SBarry Smith PetscFunctionReturn(0); 1027789d8953SBarry Smith } 1028789d8953SBarry Smith 1029789d8953SBarry Smith PetscErrorCode MatShellSetOperation_Shell(Mat mat,MatOperation op,void (*f)(void)) 1030789d8953SBarry Smith { 1031789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data;; 1032789d8953SBarry Smith 1033789d8953SBarry Smith PetscFunctionBegin; 1034789d8953SBarry Smith switch (op) { 1035789d8953SBarry Smith case MATOP_DESTROY: 1036789d8953SBarry Smith shell->ops->destroy = (PetscErrorCode (*)(Mat))f; 1037789d8953SBarry Smith break; 1038789d8953SBarry Smith case MATOP_VIEW: 1039789d8953SBarry Smith if (!mat->ops->viewnative) { 1040789d8953SBarry Smith mat->ops->viewnative = mat->ops->view; 1041789d8953SBarry Smith } 1042789d8953SBarry Smith mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 1043789d8953SBarry Smith break; 1044789d8953SBarry Smith case MATOP_COPY: 1045789d8953SBarry Smith shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f; 1046789d8953SBarry Smith break; 1047789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1048789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1049789d8953SBarry Smith case MATOP_SHIFT: 1050789d8953SBarry Smith case MATOP_SCALE: 1051789d8953SBarry Smith case MATOP_AXPY: 1052789d8953SBarry Smith case MATOP_ZERO_ROWS: 1053789d8953SBarry Smith case MATOP_ZERO_ROWS_COLUMNS: 1054789d8953SBarry Smith if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()"); 1055789d8953SBarry Smith (((void(**)(void))mat->ops)[op]) = f; 1056789d8953SBarry Smith break; 1057789d8953SBarry Smith case MATOP_GET_DIAGONAL: 1058789d8953SBarry Smith if (shell->managescalingshifts) { 1059789d8953SBarry Smith shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1060789d8953SBarry Smith mat->ops->getdiagonal = MatGetDiagonal_Shell; 1061789d8953SBarry Smith } else { 1062789d8953SBarry Smith shell->ops->getdiagonal = NULL; 1063789d8953SBarry Smith mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1064789d8953SBarry Smith } 1065789d8953SBarry Smith break; 1066789d8953SBarry Smith case MATOP_MULT: 1067789d8953SBarry Smith if (shell->managescalingshifts) { 1068789d8953SBarry Smith shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1069789d8953SBarry Smith mat->ops->mult = MatMult_Shell; 1070789d8953SBarry Smith } else { 1071789d8953SBarry Smith shell->ops->mult = NULL; 1072789d8953SBarry Smith mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1073789d8953SBarry Smith } 1074789d8953SBarry Smith break; 1075789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 1076789d8953SBarry Smith if (shell->managescalingshifts) { 1077789d8953SBarry Smith shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1078789d8953SBarry Smith mat->ops->multtranspose = MatMultTranspose_Shell; 1079789d8953SBarry Smith } else { 1080789d8953SBarry Smith shell->ops->multtranspose = NULL; 1081789d8953SBarry Smith mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1082789d8953SBarry Smith } 1083789d8953SBarry Smith break; 1084789d8953SBarry Smith default: 1085789d8953SBarry Smith (((void(**)(void))mat->ops)[op]) = f; 1086789d8953SBarry Smith break; 1087789d8953SBarry Smith } 1088789d8953SBarry Smith PetscFunctionReturn(0); 1089789d8953SBarry Smith } 1090789d8953SBarry Smith 1091789d8953SBarry Smith PetscErrorCode MatShellGetOperation_Shell(Mat mat,MatOperation op,void(**f)(void)) 1092789d8953SBarry Smith { 1093789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1094789d8953SBarry Smith 1095789d8953SBarry Smith PetscFunctionBegin; 1096789d8953SBarry Smith switch (op) { 1097789d8953SBarry Smith case MATOP_DESTROY: 1098789d8953SBarry Smith *f = (void (*)(void))shell->ops->destroy; 1099789d8953SBarry Smith break; 1100789d8953SBarry Smith case MATOP_VIEW: 1101789d8953SBarry Smith *f = (void (*)(void))mat->ops->view; 1102789d8953SBarry Smith break; 1103789d8953SBarry Smith case MATOP_COPY: 1104789d8953SBarry Smith *f = (void (*)(void))shell->ops->copy; 1105789d8953SBarry Smith break; 1106789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1107789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1108789d8953SBarry Smith case MATOP_SHIFT: 1109789d8953SBarry Smith case MATOP_SCALE: 1110789d8953SBarry Smith case MATOP_AXPY: 1111789d8953SBarry Smith case MATOP_ZERO_ROWS: 1112789d8953SBarry Smith case MATOP_ZERO_ROWS_COLUMNS: 1113789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1114789d8953SBarry Smith break; 1115789d8953SBarry Smith case MATOP_GET_DIAGONAL: 1116789d8953SBarry Smith if (shell->ops->getdiagonal) 1117789d8953SBarry Smith *f = (void (*)(void))shell->ops->getdiagonal; 1118789d8953SBarry Smith else 1119789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1120789d8953SBarry Smith break; 1121789d8953SBarry Smith case MATOP_MULT: 1122789d8953SBarry Smith if (shell->ops->mult) 1123789d8953SBarry Smith *f = (void (*)(void))shell->ops->mult; 1124789d8953SBarry Smith else 1125789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1126789d8953SBarry Smith break; 1127789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 1128789d8953SBarry Smith if (shell->ops->multtranspose) 1129789d8953SBarry Smith *f = (void (*)(void))shell->ops->multtranspose; 1130789d8953SBarry Smith else 1131789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1132789d8953SBarry Smith break; 1133789d8953SBarry Smith default: 1134789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1135789d8953SBarry Smith } 1136789d8953SBarry Smith PetscFunctionReturn(0); 1137789d8953SBarry Smith } 1138789d8953SBarry Smith 11390bad9183SKris Buschelman /*MC 1140fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 11410bad9183SKris Buschelman 11420bad9183SKris Buschelman Level: advanced 11430bad9183SKris Buschelman 11440c0fd78eSBarry Smith .seealso: MatCreateShell() 11450bad9183SKris Buschelman M*/ 11460bad9183SKris Buschelman 11478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 1148273d9f13SBarry Smith { 1149273d9f13SBarry Smith Mat_Shell *b; 1150dfbe8321SBarry Smith PetscErrorCode ierr; 1151273d9f13SBarry Smith 1152273d9f13SBarry Smith PetscFunctionBegin; 1153273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1154273d9f13SBarry Smith 1155b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1156273d9f13SBarry Smith A->data = (void*)b; 1157273d9f13SBarry Smith 115826283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 115926283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1160273d9f13SBarry Smith 1161273d9f13SBarry Smith b->ctx = 0; 1162ef66eb69SBarry Smith b->vshift = 0.0; 1163ef66eb69SBarry Smith b->vscale = 1.0; 11640c0fd78eSBarry Smith b->managescalingshifts = PETSC_TRUE; 1165273d9f13SBarry Smith A->assembled = PETSC_TRUE; 1166210f0121SBarry Smith A->preallocated = PETSC_FALSE; 11672205254eSKarl Rupp 1168789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetContext_C",MatShellGetContext_Shell);CHKERRQ(ierr); 1169789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetContext_C",MatShellSetContext_Shell);CHKERRQ(ierr); 1170789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetManageScalingShifts_C",MatShellSetManageScalingShifts_Shell);CHKERRQ(ierr); 1171789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetOperation_C",MatShellSetOperation_Shell);CHKERRQ(ierr); 1172789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetOperation_C",MatShellGetOperation_Shell);CHKERRQ(ierr); 117317667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 1174273d9f13SBarry Smith PetscFunctionReturn(0); 1175273d9f13SBarry Smith } 1176e51e0e81SBarry Smith 11774b828684SBarry Smith /*@C 1178052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 1179ff756334SLois Curfman McInnes private data storage format. 1180e51e0e81SBarry Smith 1181d083f849SBarry Smith Collective 1182c7fcc2eaSBarry Smith 1183e51e0e81SBarry Smith Input Parameters: 1184c7fcc2eaSBarry Smith + comm - MPI communicator 1185c7fcc2eaSBarry Smith . m - number of local rows (must be given) 1186c7fcc2eaSBarry Smith . n - number of local columns (must be given) 1187c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 1188c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 1189c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 1190e51e0e81SBarry Smith 1191ff756334SLois Curfman McInnes Output Parameter: 119244cd7ae7SLois Curfman McInnes . A - the matrix 1193e51e0e81SBarry Smith 1194ff2fd236SBarry Smith Level: advanced 1195ff2fd236SBarry Smith 1196f39d1f56SLois Curfman McInnes Usage: 1197*5bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1198f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 1199c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1200f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 1201f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 1202f39d1f56SLois Curfman McInnes 1203ff756334SLois Curfman McInnes Notes: 1204ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 1205ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 1206ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 1207e51e0e81SBarry Smith 120895452b02SPatrick Sanan Fortran Notes: 120995452b02SPatrick Sanan To use this from Fortran with a ctx you must write an interface definition for this 1210daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 1211daf670e6SBarry Smith in as the ctx argument. 1212f38a8d56SBarry Smith 1213f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 1214f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 1215645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 1216645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 1217645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 1218645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 1219645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 1220645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 1221645985a0SLois Curfman McInnes For example, 1222f39d1f56SLois Curfman McInnes 1223f39d1f56SLois Curfman McInnes $ 1224f39d1f56SLois Curfman McInnes $ Vec x, y 1225*5bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1226645985a0SLois Curfman McInnes $ Mat A 1227f39d1f56SLois Curfman McInnes $ 1228c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 1229c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 1230f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 1231c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 1232c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 1233c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1234645985a0SLois Curfman McInnes $ MatMult(A,x,y); 123545960306SStefano Zampini $ MatDestroy(&A); 123645960306SStefano Zampini $ VecDestroy(&y); 123745960306SStefano Zampini $ VecDestroy(&x); 1238645985a0SLois Curfman McInnes $ 1239e51e0e81SBarry Smith 12409b53d762SBarry Smith 124145960306SStefano Zampini MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), MatScale(), MatZeroRows() and MatZeroRowsColumns() internally, so these 12429b53d762SBarry Smith operations cannot be overwritten unless MatShellSetManageScalingShifts() is called. 12439b53d762SBarry Smith 12449b53d762SBarry Smith 12450c0fd78eSBarry Smith For rectangular matrices do all the scalings and shifts make sense? 12460c0fd78eSBarry Smith 124795452b02SPatrick Sanan Developers Notes: 124895452b02SPatrick Sanan Regarding shifting and scaling. The general form is 12490c0fd78eSBarry Smith 12500c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 12510c0fd78eSBarry Smith 12520c0fd78eSBarry Smith The order you apply the operations is important. For example if you have a dshift then 12530c0fd78eSBarry Smith apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift 12540c0fd78eSBarry Smith you get s*vscale*A + diag(shift) 12550c0fd78eSBarry Smith 12560c0fd78eSBarry Smith A is the user provided function. 12570c0fd78eSBarry Smith 1258ad2f5c3fSBarry Smith KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for 1259ad2f5c3fSBarry Smith for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger 1260ad2f5c3fSBarry Smith an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat); 1261ad2f5c3fSBarry Smith each time the MATSHELL matrix has changed. 1262ad2f5c3fSBarry Smith 1263ad2f5c3fSBarry Smith Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided 1264ad2f5c3fSBarry Smith with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale(). 1265ad2f5c3fSBarry Smith 12660c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts() 1267e51e0e81SBarry Smith @*/ 12687087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 1269e51e0e81SBarry Smith { 1270dfbe8321SBarry Smith PetscErrorCode ierr; 1271ed3cc1f0SBarry Smith 12723a40ed3dSBarry Smith PetscFunctionBegin; 1273f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1274f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1275273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 1276273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 12770fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 1278273d9f13SBarry Smith PetscFunctionReturn(0); 1279c7fcc2eaSBarry Smith } 1280c7fcc2eaSBarry Smith 1281789d8953SBarry Smith 1282c6866cfdSSatish Balay /*@ 1283273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 1284c7fcc2eaSBarry Smith 12853f9fe445SBarry Smith Logically Collective on Mat 1286c7fcc2eaSBarry Smith 1287273d9f13SBarry Smith Input Parameters: 1288273d9f13SBarry Smith + mat - the shell matrix 1289273d9f13SBarry Smith - ctx - the context 1290273d9f13SBarry Smith 1291273d9f13SBarry Smith Level: advanced 1292273d9f13SBarry Smith 129395452b02SPatrick Sanan Fortran Notes: 129495452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 1295daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 1296273d9f13SBarry Smith 1297273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 12980bc0a719SBarry Smith @*/ 12997087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 1300273d9f13SBarry Smith { 1301dfbe8321SBarry Smith PetscErrorCode ierr; 1302273d9f13SBarry Smith 1303273d9f13SBarry Smith PetscFunctionBegin; 13040700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1305789d8953SBarry Smith ierr = PetscUseMethod(mat,"MatShellSetContext_C",(Mat,void*),(mat,ctx));CHKERRQ(ierr); 13063a40ed3dSBarry Smith PetscFunctionReturn(0); 1307e51e0e81SBarry Smith } 1308e51e0e81SBarry Smith 13090c0fd78eSBarry Smith /*@ 13100c0fd78eSBarry Smith MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately 13110c0fd78eSBarry Smith after MatCreateShell() 13120c0fd78eSBarry Smith 13130c0fd78eSBarry Smith Logically Collective on Mat 13140c0fd78eSBarry Smith 13150c0fd78eSBarry Smith Input Parameter: 13160c0fd78eSBarry Smith . mat - the shell matrix 13170c0fd78eSBarry Smith 13180c0fd78eSBarry Smith Level: advanced 13190c0fd78eSBarry Smith 13200c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation() 13210c0fd78eSBarry Smith @*/ 13220c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A) 13230c0fd78eSBarry Smith { 13240c0fd78eSBarry Smith PetscErrorCode ierr; 13250c0fd78eSBarry Smith 13260c0fd78eSBarry Smith PetscFunctionBegin; 13270c0fd78eSBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1328789d8953SBarry Smith ierr = PetscUseMethod(A,"MatShellSetManageScalingShifts_C",(Mat),(A));CHKERRQ(ierr); 13290c0fd78eSBarry Smith PetscFunctionReturn(0); 13300c0fd78eSBarry Smith } 13310c0fd78eSBarry Smith 1332c16cb8f2SBarry Smith /*@C 1333f3b1f45cSBarry Smith MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function. 1334f3b1f45cSBarry Smith 1335f3b1f45cSBarry Smith Logically Collective on Mat 1336f3b1f45cSBarry Smith 1337f3b1f45cSBarry Smith Input Parameters: 1338f3b1f45cSBarry Smith + mat - the shell matrix 1339f3b1f45cSBarry Smith . f - the function 1340f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1341f3b1f45cSBarry Smith - ctx - an optional context for the function 1342f3b1f45cSBarry Smith 1343f3b1f45cSBarry Smith Output Parameter: 1344f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1345f3b1f45cSBarry Smith 1346f3b1f45cSBarry Smith Options Database: 1347f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1348f3b1f45cSBarry Smith 1349f3b1f45cSBarry Smith Level: advanced 1350f3b1f45cSBarry Smith 135195452b02SPatrick Sanan Fortran Notes: 135295452b02SPatrick Sanan Not supported from Fortran 1353f3b1f45cSBarry Smith 1354f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose() 1355f3b1f45cSBarry Smith @*/ 1356f3b1f45cSBarry Smith PetscErrorCode MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1357f3b1f45cSBarry Smith { 1358f3b1f45cSBarry Smith PetscErrorCode ierr; 1359f3b1f45cSBarry Smith PetscInt m,n; 1360f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1361f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 136274e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1363f3b1f45cSBarry Smith 1364f3b1f45cSBarry Smith PetscFunctionBegin; 1365f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1366f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr); 1367f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 1368f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr); 1369f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 1370f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 1371f3b1f45cSBarry Smith 13720bacdadaSStefano Zampini ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr); 13730bacdadaSStefano Zampini ierr = MatComputeOperator(mat,MATAIJ,&Dmat);CHKERRQ(ierr); 1374f3b1f45cSBarry Smith 1375f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 1376f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 1377f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 1378f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 1379f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1380f3b1f45cSBarry Smith flag = PETSC_FALSE; 1381f3b1f45cSBarry Smith if (v) { 1382fc7aafd1SBarry 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); 1383f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 1384f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 1385f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 1386f3b1f45cSBarry Smith } 1387f3b1f45cSBarry Smith } else if (v) { 1388fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr); 1389f3b1f45cSBarry Smith } 1390f3b1f45cSBarry Smith if (flg) *flg = flag; 1391f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 1392f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 1393f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 1394f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 1395f3b1f45cSBarry Smith PetscFunctionReturn(0); 1396f3b1f45cSBarry Smith } 1397f3b1f45cSBarry Smith 1398f3b1f45cSBarry Smith /*@C 1399f3b1f45cSBarry Smith MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function. 1400f3b1f45cSBarry Smith 1401f3b1f45cSBarry Smith Logically Collective on Mat 1402f3b1f45cSBarry Smith 1403f3b1f45cSBarry Smith Input Parameters: 1404f3b1f45cSBarry Smith + mat - the shell matrix 1405f3b1f45cSBarry Smith . f - the function 1406f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1407f3b1f45cSBarry Smith - ctx - an optional context for the function 1408f3b1f45cSBarry Smith 1409f3b1f45cSBarry Smith Output Parameter: 1410f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1411f3b1f45cSBarry Smith 1412f3b1f45cSBarry Smith Options Database: 1413f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1414f3b1f45cSBarry Smith 1415f3b1f45cSBarry Smith Level: advanced 1416f3b1f45cSBarry Smith 141795452b02SPatrick Sanan Fortran Notes: 141895452b02SPatrick Sanan Not supported from Fortran 1419f3b1f45cSBarry Smith 1420f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult() 1421f3b1f45cSBarry Smith @*/ 1422f3b1f45cSBarry Smith PetscErrorCode MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1423f3b1f45cSBarry Smith { 1424f3b1f45cSBarry Smith PetscErrorCode ierr; 1425f3b1f45cSBarry Smith Vec x,y,z; 1426f3b1f45cSBarry Smith PetscInt m,n,M,N; 1427f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1428f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 142974e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1430f3b1f45cSBarry Smith 1431f3b1f45cSBarry Smith PetscFunctionBegin; 1432f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1433f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr); 1434f3b1f45cSBarry Smith ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr); 1435f3b1f45cSBarry Smith ierr = VecDuplicate(y,&z);CHKERRQ(ierr); 1436f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 1437f3b1f45cSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 1438f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr); 1439f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 1440f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 14410bacdadaSStefano Zampini ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr); 1442f3b1f45cSBarry Smith ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr); 14430bacdadaSStefano Zampini ierr = MatComputeOperatorTranspose(mat,MATAIJ,&Dmat);CHKERRQ(ierr); 1444f3b1f45cSBarry Smith 1445f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 1446f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 1447f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 1448f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 1449f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1450f3b1f45cSBarry Smith flag = PETSC_FALSE; 1451f3b1f45cSBarry Smith if (v) { 1452fc7aafd1SBarry 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); 1453f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 1454f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 1455f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 1456f3b1f45cSBarry Smith } 1457f3b1f45cSBarry Smith } else if (v) { 1458fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr); 1459f3b1f45cSBarry Smith } 1460f3b1f45cSBarry Smith if (flg) *flg = flag; 1461f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 1462f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 1463f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 1464f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 1465f3b1f45cSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1466f3b1f45cSBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 1467f3b1f45cSBarry Smith ierr = VecDestroy(&z);CHKERRQ(ierr); 1468f3b1f45cSBarry Smith PetscFunctionReturn(0); 1469f3b1f45cSBarry Smith } 1470f3b1f45cSBarry Smith 1471f3b1f45cSBarry Smith /*@C 14720c0fd78eSBarry Smith MatShellSetOperation - Allows user to set a matrix operation for a shell matrix. 1473e51e0e81SBarry Smith 14743f9fe445SBarry Smith Logically Collective on Mat 1475fee21e36SBarry Smith 1476c7fcc2eaSBarry Smith Input Parameters: 1477c7fcc2eaSBarry Smith + mat - the shell matrix 1478c7fcc2eaSBarry Smith . op - the name of the operation 1479789d8953SBarry Smith - g - the function that provides the operation. 1480c7fcc2eaSBarry Smith 148115091d37SBarry Smith Level: advanced 148215091d37SBarry Smith 1483fae171e0SBarry Smith Usage: 1484e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 1485f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 1486c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 14870b627109SLois Curfman McInnes 1488a62d957aSLois Curfman McInnes Notes: 1489e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 14901c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 1491a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 14921c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 1493a62d957aSLois Curfman McInnes 1494a4d85fd7SBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 1495deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 1496deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 1497deebb3c3SLois Curfman McInnes routines, e.g., 1498a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1499a62d957aSLois Curfman McInnes 15004aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 15014aa34b0aSBarry Smith nonzero on failure. 15024aa34b0aSBarry Smith 1503a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 1504a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 1505a62d957aSLois Curfman McInnes set by MatCreateShell(). 1506a62d957aSLois Curfman McInnes 150795452b02SPatrick Sanan Fortran Notes: 150895452b02SPatrick Sanan For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 1509501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 1510501d9185SBarry Smith 15110c0fd78eSBarry Smith Use MatSetOperation() to set an operation for any matrix type 15120c0fd78eSBarry Smith 15130c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts() 1514e51e0e81SBarry Smith @*/ 1515789d8953SBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*g)(void)) 1516e51e0e81SBarry Smith { 1517976e8c5aSLisandro Dalcin PetscErrorCode ierr; 1518273d9f13SBarry Smith 15193a40ed3dSBarry Smith PetscFunctionBegin; 15200700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1521789d8953SBarry Smith ierr = PetscUseMethod(mat,"MatShellSetOperation_C",(Mat,MatOperation,void (*)(void)),(mat,op,g));CHKERRQ(ierr); 15223a40ed3dSBarry Smith PetscFunctionReturn(0); 1523e51e0e81SBarry Smith } 1524f0479e8cSBarry Smith 1525d4bb536fSBarry Smith /*@C 1526d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 1527d4bb536fSBarry Smith 1528c7fcc2eaSBarry Smith Not Collective 1529c7fcc2eaSBarry Smith 1530d4bb536fSBarry Smith Input Parameters: 1531c7fcc2eaSBarry Smith + mat - the shell matrix 1532c7fcc2eaSBarry Smith - op - the name of the operation 1533d4bb536fSBarry Smith 1534d4bb536fSBarry Smith Output Parameter: 1535789d8953SBarry Smith . g - the function that provides the operation. 1536d4bb536fSBarry Smith 153715091d37SBarry Smith Level: advanced 153815091d37SBarry Smith 1539d4bb536fSBarry Smith Notes: 1540e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 1541d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 1542d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 1543d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 1544d4bb536fSBarry Smith 1545d4bb536fSBarry Smith All user-provided functions have the same calling 1546d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 1547d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 1548d4bb536fSBarry Smith routines, e.g., 1549d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1550d4bb536fSBarry Smith 1551d4bb536fSBarry Smith Within each user-defined routine, the user should call 1552d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 1553d4bb536fSBarry Smith set by MatCreateShell(). 1554d4bb536fSBarry Smith 1555ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 1556d4bb536fSBarry Smith @*/ 1557789d8953SBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**g)(void)) 1558d4bb536fSBarry Smith { 1559976e8c5aSLisandro Dalcin PetscErrorCode ierr; 1560273d9f13SBarry Smith 15613a40ed3dSBarry Smith PetscFunctionBegin; 15620700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1563789d8953SBarry Smith ierr = PetscUseMethod(mat,"MatShellGetOperation_C",(Mat,MatOperation,void (**)(void)),(mat,op,g));CHKERRQ(ierr); 15643a40ed3dSBarry Smith PetscFunctionReturn(0); 1565d4bb536fSBarry Smith } 1566