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 2259d225801SJed Brown /*@ 226a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 227b4fd4287SBarry Smith 228c7fcc2eaSBarry Smith Not Collective 229c7fcc2eaSBarry Smith 230b4fd4287SBarry Smith Input Parameter: 231b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 232b4fd4287SBarry Smith 233b4fd4287SBarry Smith Output Parameter: 234b4fd4287SBarry Smith . ctx - the user provided context 235b4fd4287SBarry Smith 23615091d37SBarry Smith Level: advanced 23715091d37SBarry Smith 23895452b02SPatrick Sanan Fortran Notes: 23995452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 240daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 241a62d957aSLois Curfman McInnes 242ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 2430bc0a719SBarry Smith @*/ 2448fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 245b4fd4287SBarry Smith { 246dfbe8321SBarry Smith PetscErrorCode ierr; 247ace3abfcSBarry Smith PetscBool flg; 248273d9f13SBarry Smith 2493a40ed3dSBarry Smith PetscFunctionBegin; 2500700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2514482741eSBarry Smith PetscValidPointer(ctx,2); 252251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 253940183f3SBarry Smith if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 254ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix"); 2553a40ed3dSBarry Smith PetscFunctionReturn(0); 256b4fd4287SBarry Smith } 257b4fd4287SBarry Smith 25845960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Local_Shell(Mat mat,PetscInt nr,PetscInt rows[],PetscInt nc,PetscInt cols[],PetscScalar diag,PetscBool rc) 25945960306SStefano Zampini { 26045960306SStefano Zampini PetscErrorCode ierr; 26145960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)mat->data; 26245960306SStefano Zampini Vec x = NULL,b = NULL; 26345960306SStefano Zampini IS is1, is2; 26445960306SStefano Zampini const PetscInt *ridxs; 26545960306SStefano Zampini PetscInt *idxs,*gidxs; 26645960306SStefano Zampini PetscInt cum,rst,cst,i; 26745960306SStefano Zampini 26845960306SStefano Zampini PetscFunctionBegin; 26945960306SStefano Zampini if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation"); 27045960306SStefano Zampini if (!shell->zvals) { 27145960306SStefano Zampini ierr = MatCreateVecs(mat,NULL,&shell->zvals);CHKERRQ(ierr); 27245960306SStefano Zampini } 27345960306SStefano Zampini if (!shell->zvals_w) { 27445960306SStefano Zampini ierr = VecDuplicate(shell->zvals,&shell->zvals_w);CHKERRQ(ierr); 27545960306SStefano Zampini } 27645960306SStefano Zampini ierr = MatGetOwnershipRange(mat,&rst,NULL);CHKERRQ(ierr); 27745960306SStefano Zampini ierr = MatGetOwnershipRangeColumn(mat,&cst,NULL);CHKERRQ(ierr); 27845960306SStefano Zampini 27945960306SStefano Zampini /* Expand/create index set of zeroed rows */ 28045960306SStefano Zampini ierr = PetscMalloc1(nr,&idxs);CHKERRQ(ierr); 28145960306SStefano Zampini for (i = 0; i < nr; i++) idxs[i] = rows[i] + rst; 28245960306SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,nr,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr); 28345960306SStefano Zampini ierr = ISSort(is1);CHKERRQ(ierr); 28445960306SStefano Zampini ierr = VecISSet(shell->zvals,is1,diag);CHKERRQ(ierr); 28545960306SStefano Zampini if (shell->zrows) { 28645960306SStefano Zampini ierr = ISSum(shell->zrows,is1,&is2);CHKERRQ(ierr); 28745960306SStefano Zampini ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr); 28845960306SStefano Zampini ierr = ISDestroy(&is1);CHKERRQ(ierr); 28945960306SStefano Zampini shell->zrows = is2; 29045960306SStefano Zampini } else shell->zrows = is1; 29145960306SStefano Zampini 29245960306SStefano Zampini /* Create scatters for diagonal values communications */ 29345960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr); 29445960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr); 29545960306SStefano Zampini 29645960306SStefano Zampini /* row scatter: from/to left vector */ 29745960306SStefano Zampini ierr = MatCreateVecs(mat,&x,&b);CHKERRQ(ierr); 29845960306SStefano Zampini ierr = VecScatterCreate(b,shell->zrows,shell->zvals_w,shell->zrows,&shell->zvals_sct_r);CHKERRQ(ierr); 29945960306SStefano Zampini 30045960306SStefano Zampini /* col scatter: from right vector to left vector */ 30145960306SStefano Zampini ierr = ISGetIndices(shell->zrows,&ridxs);CHKERRQ(ierr); 30245960306SStefano Zampini ierr = ISGetLocalSize(shell->zrows,&nr);CHKERRQ(ierr); 30345960306SStefano Zampini ierr = PetscMalloc1(nr,&gidxs);CHKERRQ(ierr); 30445960306SStefano Zampini for (i = 0, cum = 0; i < nr; i++) { 30545960306SStefano Zampini if (ridxs[i] >= mat->cmap->N) continue; 30645960306SStefano Zampini gidxs[cum] = ridxs[i]; 30745960306SStefano Zampini cum++; 30845960306SStefano Zampini } 30945960306SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)mat),cum,gidxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr); 31045960306SStefano Zampini ierr = VecScatterCreate(x,is1,shell->zvals_w,is1,&shell->zvals_sct_c);CHKERRQ(ierr); 31145960306SStefano Zampini ierr = ISDestroy(&is1);CHKERRQ(ierr); 31245960306SStefano Zampini ierr = VecDestroy(&x);CHKERRQ(ierr); 31345960306SStefano Zampini ierr = VecDestroy(&b);CHKERRQ(ierr); 31445960306SStefano Zampini 31545960306SStefano Zampini /* Expand/create index set of zeroed columns */ 31645960306SStefano Zampini if (rc) { 31745960306SStefano Zampini ierr = PetscMalloc1(nc,&idxs);CHKERRQ(ierr); 31845960306SStefano Zampini for (i = 0; i < nc; i++) idxs[i] = cols[i] + cst; 31945960306SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,nc,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr); 32045960306SStefano Zampini ierr = ISSort(is1);CHKERRQ(ierr); 32145960306SStefano Zampini if (shell->zcols) { 32245960306SStefano Zampini ierr = ISSum(shell->zcols,is1,&is2);CHKERRQ(ierr); 32345960306SStefano Zampini ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr); 32445960306SStefano Zampini ierr = ISDestroy(&is1);CHKERRQ(ierr); 32545960306SStefano Zampini shell->zcols = is2; 32645960306SStefano Zampini } else shell->zcols = is1; 32745960306SStefano Zampini } 32845960306SStefano Zampini PetscFunctionReturn(0); 32945960306SStefano Zampini } 33045960306SStefano Zampini 33145960306SStefano Zampini static PetscErrorCode MatZeroRows_Shell(Mat mat,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 33245960306SStefano Zampini { 33345960306SStefano Zampini PetscInt nr, *lrows; 33445960306SStefano Zampini PetscErrorCode ierr; 33545960306SStefano Zampini 33645960306SStefano Zampini PetscFunctionBegin; 33745960306SStefano Zampini if (x && b) { 33845960306SStefano Zampini Vec xt; 33945960306SStefano Zampini PetscScalar *vals; 34045960306SStefano Zampini PetscInt *gcols,i,st,nl,nc; 34145960306SStefano Zampini 34245960306SStefano Zampini ierr = PetscMalloc1(n,&gcols);CHKERRQ(ierr); 34345960306SStefano Zampini for (i = 0, nc = 0; i < n; i++) if (rows[i] < mat->cmap->N) gcols[nc++] = rows[i]; 34445960306SStefano Zampini 34545960306SStefano Zampini ierr = MatCreateVecs(mat,&xt,NULL);CHKERRQ(ierr); 34645960306SStefano Zampini ierr = VecCopy(x,xt);CHKERRQ(ierr); 34745960306SStefano Zampini ierr = PetscCalloc1(nc,&vals);CHKERRQ(ierr); 34845960306SStefano Zampini ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */ 34945960306SStefano Zampini ierr = PetscFree(vals);CHKERRQ(ierr); 35045960306SStefano Zampini ierr = VecAssemblyBegin(xt);CHKERRQ(ierr); 35145960306SStefano Zampini ierr = VecAssemblyEnd(xt);CHKERRQ(ierr); 35245960306SStefano Zampini ierr = VecAYPX(xt,-1.0,x);CHKERRQ(ierr); /* xt = [0, x2] */ 35345960306SStefano Zampini 35445960306SStefano Zampini ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr); 35545960306SStefano Zampini ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr); 35645960306SStefano Zampini ierr = VecGetArray(xt,&vals);CHKERRQ(ierr); 35745960306SStefano Zampini for (i = 0; i < nl; i++) { 35845960306SStefano Zampini PetscInt g = i + st; 35945960306SStefano Zampini if (g > mat->rmap->N) continue; 36045960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 36145960306SStefano Zampini ierr = VecSetValue(b,g,diag*vals[i],INSERT_VALUES);CHKERRQ(ierr); 36245960306SStefano Zampini } 36345960306SStefano Zampini ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr); 36445960306SStefano Zampini ierr = VecAssemblyBegin(b);CHKERRQ(ierr); 36545960306SStefano Zampini ierr = VecAssemblyEnd(b);CHKERRQ(ierr); /* b = [b1, x2 * diag] */ 36645960306SStefano Zampini ierr = VecDestroy(&xt);CHKERRQ(ierr); 36745960306SStefano Zampini ierr = PetscFree(gcols);CHKERRQ(ierr); 36845960306SStefano Zampini } 36945960306SStefano Zampini ierr = PetscLayoutMapLocal(mat->rmap,n,rows,&nr,&lrows,NULL);CHKERRQ(ierr); 37045960306SStefano Zampini ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,0,NULL,diag,PETSC_FALSE);CHKERRQ(ierr); 37145960306SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 37245960306SStefano Zampini PetscFunctionReturn(0); 37345960306SStefano Zampini } 37445960306SStefano Zampini 37545960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Shell(Mat mat,PetscInt n,const PetscInt rowscols[],PetscScalar diag,Vec x,Vec b) 37645960306SStefano Zampini { 37745960306SStefano Zampini PetscInt *lrows, *lcols; 37845960306SStefano Zampini PetscInt nr, nc; 37945960306SStefano Zampini PetscBool congruent; 38045960306SStefano Zampini PetscErrorCode ierr; 38145960306SStefano Zampini 38245960306SStefano Zampini PetscFunctionBegin; 38345960306SStefano Zampini if (x && b) { 38445960306SStefano Zampini Vec xt, bt; 38545960306SStefano Zampini PetscScalar *vals; 38645960306SStefano Zampini PetscInt *grows,*gcols,i,st,nl; 38745960306SStefano Zampini 38845960306SStefano Zampini ierr = PetscMalloc2(n,&grows,n,&gcols);CHKERRQ(ierr); 38945960306SStefano Zampini for (i = 0, nr = 0; i < n; i++) if (rowscols[i] < mat->rmap->N) grows[nr++] = rowscols[i]; 39045960306SStefano Zampini for (i = 0, nc = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) gcols[nc++] = rowscols[i]; 39145960306SStefano Zampini ierr = PetscCalloc1(n,&vals);CHKERRQ(ierr); 39245960306SStefano Zampini 39345960306SStefano Zampini ierr = MatCreateVecs(mat,&xt,&bt);CHKERRQ(ierr); 39445960306SStefano Zampini ierr = VecCopy(x,xt);CHKERRQ(ierr); 39545960306SStefano Zampini ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */ 39645960306SStefano Zampini ierr = VecAssemblyBegin(xt);CHKERRQ(ierr); 39745960306SStefano Zampini ierr = VecAssemblyEnd(xt);CHKERRQ(ierr); 39845960306SStefano Zampini ierr = VecAXPY(xt,-1.0,x);CHKERRQ(ierr); /* xt = [0, -x2] */ 39945960306SStefano Zampini ierr = MatMult(mat,xt,bt);CHKERRQ(ierr); /* bt = [-A12*x2,-A22*x2] */ 40045960306SStefano Zampini ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* bt = [-A12*x2,0] */ 40145960306SStefano Zampini ierr = VecAssemblyBegin(bt);CHKERRQ(ierr); 40245960306SStefano Zampini ierr = VecAssemblyEnd(bt);CHKERRQ(ierr); 40345960306SStefano Zampini ierr = VecAXPY(b,1.0,bt);CHKERRQ(ierr); /* b = [b1 - A12*x2, b2] */ 40445960306SStefano Zampini ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* b = [b1 - A12*x2, 0] */ 40545960306SStefano Zampini ierr = VecAssemblyBegin(bt);CHKERRQ(ierr); 40645960306SStefano Zampini ierr = VecAssemblyEnd(bt);CHKERRQ(ierr); 40745960306SStefano Zampini ierr = PetscFree(vals);CHKERRQ(ierr); 40845960306SStefano Zampini 40945960306SStefano Zampini ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr); 41045960306SStefano Zampini ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr); 41145960306SStefano Zampini ierr = VecGetArray(xt,&vals);CHKERRQ(ierr); 41245960306SStefano Zampini for (i = 0; i < nl; i++) { 41345960306SStefano Zampini PetscInt g = i + st; 41445960306SStefano Zampini if (g > mat->rmap->N) continue; 41545960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 41645960306SStefano Zampini ierr = VecSetValue(b,g,-diag*vals[i],INSERT_VALUES);CHKERRQ(ierr); 41745960306SStefano Zampini } 41845960306SStefano Zampini ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr); 41945960306SStefano Zampini ierr = VecAssemblyBegin(b);CHKERRQ(ierr); 42045960306SStefano Zampini ierr = VecAssemblyEnd(b);CHKERRQ(ierr); /* b = [b1 - A12*x2, x2 * diag] */ 42145960306SStefano Zampini ierr = VecDestroy(&xt);CHKERRQ(ierr); 42245960306SStefano Zampini ierr = VecDestroy(&bt);CHKERRQ(ierr); 42345960306SStefano Zampini ierr = PetscFree2(grows,gcols);CHKERRQ(ierr); 42445960306SStefano Zampini } 42545960306SStefano Zampini ierr = PetscLayoutMapLocal(mat->rmap,n,rowscols,&nr,&lrows,NULL);CHKERRQ(ierr); 42645960306SStefano Zampini ierr = MatHasCongruentLayouts(mat,&congruent);CHKERRQ(ierr); 42745960306SStefano Zampini if (congruent) { 42845960306SStefano Zampini nc = nr; 42945960306SStefano Zampini lcols = lrows; 43045960306SStefano Zampini } else { /* MatZeroRowsColumns implicitly assumes the rowscols indices are for a square matrix, here we handle a more general case */ 43145960306SStefano Zampini PetscInt i,nt,*t; 43245960306SStefano Zampini 43345960306SStefano Zampini ierr = PetscMalloc1(n,&t);CHKERRQ(ierr); 43445960306SStefano Zampini for (i = 0, nt = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) t[nt++] = rowscols[i]; 43545960306SStefano Zampini ierr = PetscLayoutMapLocal(mat->cmap,nt,t,&nc,&lcols,NULL);CHKERRQ(ierr); 43645960306SStefano Zampini ierr = PetscFree(t);CHKERRQ(ierr); 43745960306SStefano Zampini } 43845960306SStefano Zampini ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,nc,lcols,diag,PETSC_TRUE);CHKERRQ(ierr); 43945960306SStefano Zampini if (!congruent) { 44045960306SStefano Zampini ierr = PetscFree(lcols);CHKERRQ(ierr); 44145960306SStefano Zampini } 44245960306SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 44345960306SStefano Zampini PetscFunctionReturn(0); 44445960306SStefano Zampini } 44545960306SStefano Zampini 446dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 447e51e0e81SBarry Smith { 448dfbe8321SBarry Smith PetscErrorCode ierr; 449bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 450ed3cc1f0SBarry Smith 4513a40ed3dSBarry Smith PetscFunctionBegin; 45228f357e3SAlex Fikl if (shell->ops->destroy) { 45328f357e3SAlex Fikl ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr); 454bf0cc555SLisandro Dalcin } 4550c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 4560c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 4570c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 4588fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 4598fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 4605edf6516SJed Brown ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr); 4615edf6516SJed Brown ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr); 4629f137db4SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 46345960306SStefano Zampini ierr = VecDestroy(&shell->zvals_w);CHKERRQ(ierr); 46445960306SStefano Zampini ierr = VecDestroy(&shell->zvals);CHKERRQ(ierr); 46545960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr); 46645960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr); 46745960306SStefano Zampini ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr); 46845960306SStefano Zampini ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr); 469bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 4703a40ed3dSBarry Smith PetscFunctionReturn(0); 471e51e0e81SBarry Smith } 472e51e0e81SBarry Smith 4737fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str) 4747fabda3fSAlex Fikl { 47528f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data; 4767fabda3fSAlex Fikl PetscErrorCode ierr; 477cb8c8a77SBarry Smith PetscBool matflg; 4787fabda3fSAlex Fikl 4797fabda3fSAlex Fikl PetscFunctionBegin; 48028f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr); 481cb8c8a77SBarry Smith if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL"); 4827fabda3fSAlex Fikl 483cb8c8a77SBarry Smith ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 484cb8c8a77SBarry Smith ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr); 4857fabda3fSAlex Fikl 486cb8c8a77SBarry Smith if (shellA->ops->copy) { 48728f357e3SAlex Fikl ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr); 488cb8c8a77SBarry Smith } 4897fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 4907fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 4910c0fd78eSBarry Smith if (shellA->dshift) { 4920c0fd78eSBarry Smith if (!shellB->dshift) { 4930c0fd78eSBarry Smith ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr); 4947fabda3fSAlex Fikl } 4950c0fd78eSBarry Smith ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr); 4967fabda3fSAlex Fikl } else { 4979f137db4SBarry Smith ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr); 4987fabda3fSAlex Fikl } 4990c0fd78eSBarry Smith if (shellA->left) { 5000c0fd78eSBarry Smith if (!shellB->left) { 5010c0fd78eSBarry Smith ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr); 5027fabda3fSAlex Fikl } 5030c0fd78eSBarry Smith ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr); 5047fabda3fSAlex Fikl } else { 5059f137db4SBarry Smith ierr = VecDestroy(&shellB->left);CHKERRQ(ierr); 5067fabda3fSAlex Fikl } 5070c0fd78eSBarry Smith if (shellA->right) { 5080c0fd78eSBarry Smith if (!shellB->right) { 5090c0fd78eSBarry Smith ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr); 5107fabda3fSAlex Fikl } 5110c0fd78eSBarry Smith ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr); 5127fabda3fSAlex Fikl } else { 5139f137db4SBarry Smith ierr = VecDestroy(&shellB->right);CHKERRQ(ierr); 5147fabda3fSAlex Fikl } 51540e381c3SBarry Smith ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr); 51640e381c3SBarry Smith if (shellA->axpy) { 51740e381c3SBarry Smith ierr = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr); 51840e381c3SBarry Smith shellB->axpy = shellA->axpy; 51940e381c3SBarry Smith shellB->axpy_vscale = shellA->axpy_vscale; 52040e381c3SBarry Smith } 52145960306SStefano Zampini if (shellA->zrows) { 52245960306SStefano Zampini ierr = ISDuplicate(shellA->zrows,&shellB->zrows);CHKERRQ(ierr); 52345960306SStefano Zampini if (shellA->zcols) { 52445960306SStefano Zampini ierr = ISDuplicate(shellA->zcols,&shellB->zcols);CHKERRQ(ierr); 52545960306SStefano Zampini } 52645960306SStefano Zampini ierr = VecDuplicate(shellA->zvals,&shellB->zvals);CHKERRQ(ierr); 52745960306SStefano Zampini ierr = VecCopy(shellA->zvals,shellB->zvals);CHKERRQ(ierr); 52845960306SStefano Zampini ierr = VecDuplicate(shellA->zvals_w,&shellB->zvals_w);CHKERRQ(ierr); 52945960306SStefano Zampini ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_r);CHKERRQ(ierr); 53045960306SStefano Zampini ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_c);CHKERRQ(ierr); 53145960306SStefano Zampini shellB->zvals_sct_r = shellA->zvals_sct_r; 53245960306SStefano Zampini shellB->zvals_sct_c = shellA->zvals_sct_c; 53345960306SStefano Zampini } 5347fabda3fSAlex Fikl PetscFunctionReturn(0); 5357fabda3fSAlex Fikl } 5367fabda3fSAlex Fikl 537cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M) 538cb8c8a77SBarry Smith { 539cb8c8a77SBarry Smith PetscErrorCode ierr; 540cb8c8a77SBarry Smith void *ctx; 541cb8c8a77SBarry Smith 542cb8c8a77SBarry Smith PetscFunctionBegin; 543cb8c8a77SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 544cb8c8a77SBarry Smith ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr); 545a4b1380bSStefano Zampini if (op != MAT_DO_NOT_COPY_VALUES) { 546cb8c8a77SBarry Smith ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 547a4b1380bSStefano Zampini } 548cb8c8a77SBarry Smith PetscFunctionReturn(0); 549cb8c8a77SBarry Smith } 550cb8c8a77SBarry Smith 551dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 552ef66eb69SBarry Smith { 553ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 554dfbe8321SBarry Smith PetscErrorCode ierr; 55525578ef6SJed Brown Vec xx; 556e3f487b0SBarry Smith PetscObjectState instate,outstate; 557ef66eb69SBarry Smith 558ef66eb69SBarry Smith PetscFunctionBegin; 559976e8c5aSLisandro Dalcin if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL"); 56045960306SStefano Zampini ierr = MatShellPreZeroRight(A,x,&xx);CHKERRQ(ierr); 56145960306SStefano Zampini ierr = MatShellPreScaleRight(A,xx,&xx);CHKERRQ(ierr); 56245960306SStefano Zampini ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 56328f357e3SAlex Fikl ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr); 564e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 565e3f487b0SBarry Smith if (instate == outstate) { 566e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 567e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 568e3f487b0SBarry Smith } 5698fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 5708fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 57145960306SStefano Zampini ierr = MatShellPostZeroLeft(A,y);CHKERRQ(ierr); 5729f137db4SBarry Smith 5739f137db4SBarry Smith if (shell->axpy) { 5749f137db4SBarry Smith if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);} 5759f137db4SBarry Smith ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr); 5769f137db4SBarry Smith ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 5779f137db4SBarry Smith } 578ef66eb69SBarry Smith PetscFunctionReturn(0); 579ef66eb69SBarry Smith } 580ef66eb69SBarry Smith 5815edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 5825edf6516SJed Brown { 5835edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 5845edf6516SJed Brown PetscErrorCode ierr; 5855edf6516SJed Brown 5865edf6516SJed Brown PetscFunctionBegin; 5875edf6516SJed Brown if (y == z) { 5885edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 5895edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 590b8430d49SBarry Smith ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr); 5915edf6516SJed Brown } else { 5925edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 5935edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 5945edf6516SJed Brown } 5955edf6516SJed Brown PetscFunctionReturn(0); 5965edf6516SJed Brown } 5975edf6516SJed Brown 59818be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 59918be62a5SSatish Balay { 60018be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 60118be62a5SSatish Balay PetscErrorCode ierr; 60225578ef6SJed Brown Vec xx; 603e3f487b0SBarry Smith PetscObjectState instate,outstate; 60418be62a5SSatish Balay 60518be62a5SSatish Balay PetscFunctionBegin; 60645960306SStefano Zampini ierr = MatShellPreZeroLeft(A,x,&xx);CHKERRQ(ierr); 60745960306SStefano Zampini ierr = MatShellPreScaleLeft(A,xx,&xx);CHKERRQ(ierr); 608e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 6090c0fd78eSBarry Smith if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL"); 61028f357e3SAlex Fikl ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr); 611e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 612e3f487b0SBarry Smith if (instate == outstate) { 613e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 614e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 615e3f487b0SBarry Smith } 6168fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 6178fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 61845960306SStefano Zampini ierr = MatShellPostZeroRight(A,y);CHKERRQ(ierr); 619350ec83bSStefano Zampini 620350ec83bSStefano Zampini if (shell->axpy) { 621350ec83bSStefano Zampini if (!shell->right_work) {ierr = MatCreateVecs(A,NULL,&shell->right_work);CHKERRQ(ierr);} 622350ec83bSStefano Zampini ierr = MatMultTranspose(shell->axpy,x,shell->right_work);CHKERRQ(ierr); 623350ec83bSStefano Zampini ierr = VecAXPY(y,shell->axpy_vscale,shell->right_work);CHKERRQ(ierr); 624350ec83bSStefano Zampini } 62518be62a5SSatish Balay PetscFunctionReturn(0); 62618be62a5SSatish Balay } 62718be62a5SSatish Balay 6285edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 6295edf6516SJed Brown { 6305edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 6315edf6516SJed Brown PetscErrorCode ierr; 6325edf6516SJed Brown 6335edf6516SJed Brown PetscFunctionBegin; 6345edf6516SJed Brown if (y == z) { 6355edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 6365edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 637dac989eaS“Marek ierr = VecAXPY(z,1.0,shell->left_add_work);CHKERRQ(ierr); 6385edf6516SJed Brown } else { 6395edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 6405edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 6415edf6516SJed Brown } 6425edf6516SJed Brown PetscFunctionReturn(0); 6435edf6516SJed Brown } 6445edf6516SJed Brown 6450c0fd78eSBarry Smith /* 6460c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 6470c0fd78eSBarry Smith */ 64818be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 64918be62a5SSatish Balay { 65018be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 65118be62a5SSatish Balay PetscErrorCode ierr; 65218be62a5SSatish Balay 65318be62a5SSatish Balay PetscFunctionBegin; 6540c0fd78eSBarry Smith if (shell->ops->getdiagonal) { 65528f357e3SAlex Fikl ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr); 656305211d5SBarry 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,...)"); 65718be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 6588fe8eb27SJed Brown if (shell->dshift) { 6590c0fd78eSBarry Smith ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr); 66018be62a5SSatish Balay } 6610c0fd78eSBarry Smith ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 6628fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 6638fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 66445960306SStefano Zampini if (shell->zrows) { 66545960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 66645960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 66745960306SStefano Zampini } 66881c02519SBarry Smith if (shell->axpy) { 66981c02519SBarry Smith if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);} 67081c02519SBarry Smith ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr); 67181c02519SBarry Smith ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 67281c02519SBarry Smith } 67318be62a5SSatish Balay PetscFunctionReturn(0); 67418be62a5SSatish Balay } 67518be62a5SSatish Balay 676f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 677ef66eb69SBarry Smith { 678ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 6798fe8eb27SJed Brown PetscErrorCode ierr; 680b24ad042SBarry Smith 681ef66eb69SBarry Smith PetscFunctionBegin; 6820c0fd78eSBarry Smith if (shell->left || shell->right) { 6838fe8eb27SJed Brown if (!shell->dshift) { 6840c0fd78eSBarry Smith ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr); 6850c0fd78eSBarry Smith ierr = VecSet(shell->dshift,a);CHKERRQ(ierr); 6860c0fd78eSBarry Smith } else { 6870c0fd78eSBarry Smith if (shell->left) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 6880c0fd78eSBarry Smith if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 6899f137db4SBarry Smith ierr = VecShift(shell->dshift,a);CHKERRQ(ierr); 6900c0fd78eSBarry Smith } 6918fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 6928fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 6938fe8eb27SJed Brown } else shell->vshift += a; 69445960306SStefano Zampini if (shell->zrows) { 69545960306SStefano Zampini ierr = VecShift(shell->zvals,a);CHKERRQ(ierr); 69645960306SStefano Zampini } 697ef66eb69SBarry Smith PetscFunctionReturn(0); 698ef66eb69SBarry Smith } 699ef66eb69SBarry Smith 700b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s) 7016464bf51SAlex Fikl { 7026464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 7036464bf51SAlex Fikl PetscErrorCode ierr; 7046464bf51SAlex Fikl 7056464bf51SAlex Fikl PetscFunctionBegin; 7060c0fd78eSBarry Smith if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);} 7070c0fd78eSBarry Smith if (shell->left || shell->right) { 70892fabd1fSBarry Smith if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);} 7099f137db4SBarry Smith if (shell->left && shell->right) { 7109f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 7119f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr); 7129f137db4SBarry Smith } else if (shell->left) { 7139f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 7149f137db4SBarry Smith } else { 7159f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr); 7169f137db4SBarry Smith } 717b253ae0bS“Marek ierr = VecAXPY(shell->dshift,s,shell->right_work);CHKERRQ(ierr); 7180c0fd78eSBarry Smith } else { 719b253ae0bS“Marek ierr = VecAXPY(shell->dshift,s,D);CHKERRQ(ierr); 720b253ae0bS“Marek } 721b253ae0bS“Marek PetscFunctionReturn(0); 722b253ae0bS“Marek } 723b253ae0bS“Marek 724b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 725b253ae0bS“Marek { 72645960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 727b253ae0bS“Marek Vec d; 728b253ae0bS“Marek PetscErrorCode ierr; 729b253ae0bS“Marek 730b253ae0bS“Marek PetscFunctionBegin; 731b253ae0bS“Marek if (ins == INSERT_VALUES) { 732b253ae0bS“Marek if (!A->ops->getdiagonal) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first"); 733b253ae0bS“Marek ierr = VecDuplicate(D,&d);CHKERRQ(ierr); 734b253ae0bS“Marek ierr = MatGetDiagonal(A,d);CHKERRQ(ierr); 735b253ae0bS“Marek ierr = MatDiagonalSet_Shell_Private(A,d,-1.);CHKERRQ(ierr); 736b253ae0bS“Marek ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr); 737b253ae0bS“Marek ierr = VecDestroy(&d);CHKERRQ(ierr); 73845960306SStefano Zampini if (shell->zrows) { 73945960306SStefano Zampini ierr = VecCopy(D,shell->zvals);CHKERRQ(ierr); 74045960306SStefano Zampini } 741b253ae0bS“Marek } else { 742b253ae0bS“Marek ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr); 74345960306SStefano Zampini if (shell->zrows) { 74445960306SStefano Zampini ierr = VecAXPY(shell->zvals,1.0,D);CHKERRQ(ierr); 74545960306SStefano Zampini } 7466464bf51SAlex Fikl } 7476464bf51SAlex Fikl PetscFunctionReturn(0); 7486464bf51SAlex Fikl } 7496464bf51SAlex Fikl 750f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 751ef66eb69SBarry Smith { 752ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 7538fe8eb27SJed Brown PetscErrorCode ierr; 754b24ad042SBarry Smith 755ef66eb69SBarry Smith PetscFunctionBegin; 756f4df32b1SMatthew Knepley shell->vscale *= a; 7570c0fd78eSBarry Smith shell->vshift *= a; 7582205254eSKarl Rupp if (shell->dshift) { 7592205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 7600c0fd78eSBarry Smith } 76181c02519SBarry Smith shell->axpy_vscale *= a; 76245960306SStefano Zampini if (shell->zrows) { 76345960306SStefano Zampini ierr = VecScale(shell->zvals,a);CHKERRQ(ierr); 76445960306SStefano Zampini } 7658fe8eb27SJed Brown PetscFunctionReturn(0); 76618be62a5SSatish Balay } 7678fe8eb27SJed Brown 7688fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 7698fe8eb27SJed Brown { 7708fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 7718fe8eb27SJed Brown PetscErrorCode ierr; 7728fe8eb27SJed Brown 7738fe8eb27SJed Brown PetscFunctionBegin; 77481c02519SBarry Smith if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation"); 7758fe8eb27SJed Brown if (left) { 7760c0fd78eSBarry Smith if (!shell->left) { 7770c0fd78eSBarry Smith ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr); 7788fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 7790c0fd78eSBarry Smith } else { 7800c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 78118be62a5SSatish Balay } 78245960306SStefano Zampini if (shell->zrows) { 78345960306SStefano Zampini ierr = VecPointwiseMult(shell->zvals,shell->zvals,left);CHKERRQ(ierr); 78445960306SStefano Zampini } 785ef66eb69SBarry Smith } 7868fe8eb27SJed Brown if (right) { 7870c0fd78eSBarry Smith if (!shell->right) { 7880c0fd78eSBarry Smith ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr); 7898fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 7900c0fd78eSBarry Smith } else { 7910c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 7928fe8eb27SJed Brown } 79345960306SStefano Zampini if (shell->zrows) { 79445960306SStefano Zampini if (!shell->left_work) { 79545960306SStefano Zampini ierr = MatCreateVecs(Y,NULL,&shell->left_work);CHKERRQ(ierr); 79645960306SStefano Zampini } 79745960306SStefano Zampini ierr = VecSet(shell->zvals_w,1.0);CHKERRQ(ierr); 79845960306SStefano Zampini ierr = VecScatterBegin(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 79945960306SStefano Zampini ierr = VecScatterEnd(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 80045960306SStefano Zampini ierr = VecPointwiseMult(shell->zvals,shell->zvals,shell->zvals_w);CHKERRQ(ierr); 80145960306SStefano Zampini } 8028fe8eb27SJed Brown } 803ef66eb69SBarry Smith PetscFunctionReturn(0); 804ef66eb69SBarry Smith } 805ef66eb69SBarry Smith 806dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 807ef66eb69SBarry Smith { 808ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 8090c0fd78eSBarry Smith PetscErrorCode ierr; 810ef66eb69SBarry Smith 811ef66eb69SBarry Smith PetscFunctionBegin; 8128fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 813ef66eb69SBarry Smith shell->vshift = 0.0; 814ef66eb69SBarry Smith shell->vscale = 1.0; 8150c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 8160c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 8170c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 81881c02519SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 81945960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr); 82045960306SStefano Zampini ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr); 82145960306SStefano Zampini ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr); 82245960306SStefano Zampini ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr); 823ef66eb69SBarry Smith } 824ef66eb69SBarry Smith PetscFunctionReturn(0); 825ef66eb69SBarry Smith } 826ef66eb69SBarry Smith 8273b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 8283b49f96aSBarry Smith { 8293b49f96aSBarry Smith PetscFunctionBegin; 8303b49f96aSBarry Smith *missing = PETSC_FALSE; 8313b49f96aSBarry Smith PetscFunctionReturn(0); 8323b49f96aSBarry Smith } 8333b49f96aSBarry Smith 8349f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str) 8359f137db4SBarry Smith { 8369f137db4SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 8379f137db4SBarry Smith PetscErrorCode ierr; 8389f137db4SBarry Smith 8399f137db4SBarry Smith PetscFunctionBegin; 8409f137db4SBarry Smith ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 8419f137db4SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 8429f137db4SBarry Smith shell->axpy = X; 8439f137db4SBarry Smith shell->axpy_vscale = a; 8449f137db4SBarry Smith PetscFunctionReturn(0); 8459f137db4SBarry Smith } 8469f137db4SBarry Smith 84709dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 84820563c6bSBarry Smith 0, 84920563c6bSBarry Smith 0, 8509f137db4SBarry Smith 0, 8510c0fd78eSBarry Smith /* 4*/ MatMultAdd_Shell, 8529f137db4SBarry Smith 0, 8530c0fd78eSBarry Smith MatMultTransposeAdd_Shell, 854b951964fSBarry Smith 0, 855b951964fSBarry Smith 0, 856b951964fSBarry Smith 0, 85797304618SKris Buschelman /*10*/ 0, 858b951964fSBarry Smith 0, 859b951964fSBarry Smith 0, 860b951964fSBarry Smith 0, 861b951964fSBarry Smith 0, 86297304618SKris Buschelman /*15*/ 0, 863b951964fSBarry Smith 0, 86400849d43SBarry Smith 0, 8658fe8eb27SJed Brown MatDiagonalScale_Shell, 866b951964fSBarry Smith 0, 86797304618SKris Buschelman /*20*/ 0, 868ef66eb69SBarry Smith MatAssemblyEnd_Shell, 869b951964fSBarry Smith 0, 870b951964fSBarry Smith 0, 87145960306SStefano Zampini /*24*/ MatZeroRows_Shell, 872b951964fSBarry Smith 0, 873b951964fSBarry Smith 0, 874b951964fSBarry Smith 0, 875b951964fSBarry Smith 0, 876d519adbfSMatthew Knepley /*29*/ 0, 877b951964fSBarry Smith 0, 878273d9f13SBarry Smith 0, 879b951964fSBarry Smith 0, 880b951964fSBarry Smith 0, 881cb8c8a77SBarry Smith /*34*/ MatDuplicate_Shell, 882b951964fSBarry Smith 0, 883b951964fSBarry Smith 0, 88409dc0095SBarry Smith 0, 88509dc0095SBarry Smith 0, 8869f137db4SBarry Smith /*39*/ MatAXPY_Shell, 88709dc0095SBarry Smith 0, 88809dc0095SBarry Smith 0, 88909dc0095SBarry Smith 0, 890cb8c8a77SBarry Smith MatCopy_Shell, 891d519adbfSMatthew Knepley /*44*/ 0, 892ef66eb69SBarry Smith MatScale_Shell, 893ef66eb69SBarry Smith MatShift_Shell, 8946464bf51SAlex Fikl MatDiagonalSet_Shell, 89545960306SStefano Zampini MatZeroRowsColumns_Shell, 896f73d5cc4SBarry Smith /*49*/ 0, 89709dc0095SBarry Smith 0, 89809dc0095SBarry Smith 0, 89909dc0095SBarry Smith 0, 90009dc0095SBarry Smith 0, 901d519adbfSMatthew Knepley /*54*/ 0, 90209dc0095SBarry Smith 0, 90309dc0095SBarry Smith 0, 90409dc0095SBarry Smith 0, 90509dc0095SBarry Smith 0, 906d519adbfSMatthew Knepley /*59*/ 0, 907b9b97703SBarry Smith MatDestroy_Shell, 90809dc0095SBarry Smith 0, 909251fad3fSStefano Zampini MatConvertFrom_Shell, 910273d9f13SBarry Smith 0, 911d519adbfSMatthew Knepley /*64*/ 0, 912273d9f13SBarry Smith 0, 913273d9f13SBarry Smith 0, 914273d9f13SBarry Smith 0, 915273d9f13SBarry Smith 0, 916d519adbfSMatthew Knepley /*69*/ 0, 917273d9f13SBarry Smith 0, 918c87e5d42SMatthew Knepley MatConvert_Shell, 919273d9f13SBarry Smith 0, 92097304618SKris Buschelman 0, 921d519adbfSMatthew Knepley /*74*/ 0, 92297304618SKris Buschelman 0, 92397304618SKris Buschelman 0, 92497304618SKris Buschelman 0, 92597304618SKris Buschelman 0, 926d519adbfSMatthew Knepley /*79*/ 0, 92797304618SKris Buschelman 0, 92897304618SKris Buschelman 0, 92997304618SKris Buschelman 0, 930865e5f61SKris Buschelman 0, 931d519adbfSMatthew Knepley /*84*/ 0, 932865e5f61SKris Buschelman 0, 933865e5f61SKris Buschelman 0, 934865e5f61SKris Buschelman 0, 935865e5f61SKris Buschelman 0, 936d519adbfSMatthew Knepley /*89*/ 0, 937865e5f61SKris Buschelman 0, 938865e5f61SKris Buschelman 0, 939865e5f61SKris Buschelman 0, 940865e5f61SKris Buschelman 0, 941d519adbfSMatthew Knepley /*94*/ 0, 942865e5f61SKris Buschelman 0, 943865e5f61SKris Buschelman 0, 9443964eb88SJed Brown 0, 9453964eb88SJed Brown 0, 9463964eb88SJed Brown /*99*/ 0, 9473964eb88SJed Brown 0, 9483964eb88SJed Brown 0, 9493964eb88SJed Brown 0, 9503964eb88SJed Brown 0, 9513964eb88SJed Brown /*104*/ 0, 9523964eb88SJed Brown 0, 9533964eb88SJed Brown 0, 9543964eb88SJed Brown 0, 9553964eb88SJed Brown 0, 9563964eb88SJed Brown /*109*/ 0, 9573964eb88SJed Brown 0, 9583964eb88SJed Brown 0, 9593964eb88SJed Brown 0, 9603b49f96aSBarry Smith MatMissingDiagonal_Shell, 9613964eb88SJed Brown /*114*/ 0, 9623964eb88SJed Brown 0, 9633964eb88SJed Brown 0, 9643964eb88SJed Brown 0, 9653964eb88SJed Brown 0, 9663964eb88SJed Brown /*119*/ 0, 9673964eb88SJed Brown 0, 9683964eb88SJed Brown 0, 9693964eb88SJed Brown 0, 9703964eb88SJed Brown 0, 9713964eb88SJed Brown /*124*/ 0, 9723964eb88SJed Brown 0, 9733964eb88SJed Brown 0, 9743964eb88SJed Brown 0, 9753964eb88SJed Brown 0, 9763964eb88SJed Brown /*129*/ 0, 9773964eb88SJed Brown 0, 9783964eb88SJed Brown 0, 9793964eb88SJed Brown 0, 9803964eb88SJed Brown 0, 9813964eb88SJed Brown /*134*/ 0, 9823964eb88SJed Brown 0, 9833964eb88SJed Brown 0, 9843964eb88SJed Brown 0, 9853964eb88SJed Brown 0, 9863964eb88SJed Brown /*139*/ 0, 987f9426fe0SMark Adams 0, 9883964eb88SJed Brown 0 9893964eb88SJed Brown }; 990273d9f13SBarry Smith 9910bad9183SKris Buschelman /*MC 992fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 9930bad9183SKris Buschelman 9940bad9183SKris Buschelman Level: advanced 9950bad9183SKris Buschelman 9960c0fd78eSBarry Smith .seealso: MatCreateShell() 9970bad9183SKris Buschelman M*/ 9980bad9183SKris Buschelman 9998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 1000273d9f13SBarry Smith { 1001273d9f13SBarry Smith Mat_Shell *b; 1002dfbe8321SBarry Smith PetscErrorCode ierr; 1003273d9f13SBarry Smith 1004273d9f13SBarry Smith PetscFunctionBegin; 1005273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1006273d9f13SBarry Smith 1007b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1008273d9f13SBarry Smith A->data = (void*)b; 1009273d9f13SBarry Smith 101026283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 101126283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1012273d9f13SBarry Smith 1013273d9f13SBarry Smith b->ctx = 0; 1014ef66eb69SBarry Smith b->vshift = 0.0; 1015ef66eb69SBarry Smith b->vscale = 1.0; 10160c0fd78eSBarry Smith b->managescalingshifts = PETSC_TRUE; 1017273d9f13SBarry Smith A->assembled = PETSC_TRUE; 1018210f0121SBarry Smith A->preallocated = PETSC_FALSE; 10192205254eSKarl Rupp 102017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 1021273d9f13SBarry Smith PetscFunctionReturn(0); 1022273d9f13SBarry Smith } 1023e51e0e81SBarry Smith 10244b828684SBarry Smith /*@C 1025052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 1026ff756334SLois Curfman McInnes private data storage format. 1027e51e0e81SBarry Smith 1028*d083f849SBarry Smith Collective 1029c7fcc2eaSBarry Smith 1030e51e0e81SBarry Smith Input Parameters: 1031c7fcc2eaSBarry Smith + comm - MPI communicator 1032c7fcc2eaSBarry Smith . m - number of local rows (must be given) 1033c7fcc2eaSBarry Smith . n - number of local columns (must be given) 1034c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 1035c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 1036c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 1037e51e0e81SBarry Smith 1038ff756334SLois Curfman McInnes Output Parameter: 103944cd7ae7SLois Curfman McInnes . A - the matrix 1040e51e0e81SBarry Smith 1041ff2fd236SBarry Smith Level: advanced 1042ff2fd236SBarry Smith 1043f39d1f56SLois Curfman McInnes Usage: 10447b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 1045f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 1046c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1047f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 1048f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 1049f39d1f56SLois Curfman McInnes 1050ff756334SLois Curfman McInnes Notes: 1051ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 1052ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 1053ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 1054e51e0e81SBarry Smith 105595452b02SPatrick Sanan Fortran Notes: 105695452b02SPatrick Sanan To use this from Fortran with a ctx you must write an interface definition for this 1057daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 1058daf670e6SBarry Smith in as the ctx argument. 1059f38a8d56SBarry Smith 1060f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 1061f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 1062645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 1063645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 1064645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 1065645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 1066645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 1067645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 1068645985a0SLois Curfman McInnes For example, 1069f39d1f56SLois Curfman McInnes 1070f39d1f56SLois Curfman McInnes $ 1071f39d1f56SLois Curfman McInnes $ Vec x, y 10727b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 1073645985a0SLois Curfman McInnes $ Mat A 1074f39d1f56SLois Curfman McInnes $ 1075c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 1076c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 1077f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 1078c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 1079c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 1080c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1081645985a0SLois Curfman McInnes $ MatMult(A,x,y); 108245960306SStefano Zampini $ MatDestroy(&A); 108345960306SStefano Zampini $ VecDestroy(&y); 108445960306SStefano Zampini $ VecDestroy(&x); 1085645985a0SLois Curfman McInnes $ 1086e51e0e81SBarry Smith 10879b53d762SBarry Smith 108845960306SStefano Zampini MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), MatScale(), MatZeroRows() and MatZeroRowsColumns() internally, so these 10899b53d762SBarry Smith operations cannot be overwritten unless MatShellSetManageScalingShifts() is called. 10909b53d762SBarry Smith 10919b53d762SBarry Smith 10920c0fd78eSBarry Smith For rectangular matrices do all the scalings and shifts make sense? 10930c0fd78eSBarry Smith 109495452b02SPatrick Sanan Developers Notes: 109595452b02SPatrick Sanan Regarding shifting and scaling. The general form is 10960c0fd78eSBarry Smith 10970c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 10980c0fd78eSBarry Smith 10990c0fd78eSBarry Smith The order you apply the operations is important. For example if you have a dshift then 11000c0fd78eSBarry Smith apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift 11010c0fd78eSBarry Smith you get s*vscale*A + diag(shift) 11020c0fd78eSBarry Smith 11030c0fd78eSBarry Smith A is the user provided function. 11040c0fd78eSBarry Smith 1105ad2f5c3fSBarry Smith KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for 1106ad2f5c3fSBarry Smith for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger 1107ad2f5c3fSBarry Smith an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat); 1108ad2f5c3fSBarry Smith each time the MATSHELL matrix has changed. 1109ad2f5c3fSBarry Smith 1110ad2f5c3fSBarry Smith Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided 1111ad2f5c3fSBarry Smith with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale(). 1112ad2f5c3fSBarry Smith 11130c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts() 1114e51e0e81SBarry Smith @*/ 11157087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 1116e51e0e81SBarry Smith { 1117dfbe8321SBarry Smith PetscErrorCode ierr; 1118ed3cc1f0SBarry Smith 11193a40ed3dSBarry Smith PetscFunctionBegin; 1120f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1121f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1122273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 1123273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 11240fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 1125273d9f13SBarry Smith PetscFunctionReturn(0); 1126c7fcc2eaSBarry Smith } 1127c7fcc2eaSBarry Smith 1128c6866cfdSSatish Balay /*@ 1129273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 1130c7fcc2eaSBarry Smith 11313f9fe445SBarry Smith Logically Collective on Mat 1132c7fcc2eaSBarry Smith 1133273d9f13SBarry Smith Input Parameters: 1134273d9f13SBarry Smith + mat - the shell matrix 1135273d9f13SBarry Smith - ctx - the context 1136273d9f13SBarry Smith 1137273d9f13SBarry Smith Level: advanced 1138273d9f13SBarry Smith 113995452b02SPatrick Sanan Fortran Notes: 114095452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 1141daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 1142273d9f13SBarry Smith 1143273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 11440bc0a719SBarry Smith @*/ 11457087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 1146273d9f13SBarry Smith { 1147273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1148dfbe8321SBarry Smith PetscErrorCode ierr; 1149ace3abfcSBarry Smith PetscBool flg; 1150273d9f13SBarry Smith 1151273d9f13SBarry Smith PetscFunctionBegin; 11520700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1153251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 1154273d9f13SBarry Smith if (flg) { 1155273d9f13SBarry Smith shell->ctx = ctx; 1156ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix"); 11573a40ed3dSBarry Smith PetscFunctionReturn(0); 1158e51e0e81SBarry Smith } 1159e51e0e81SBarry Smith 11600c0fd78eSBarry Smith /*@ 11610c0fd78eSBarry Smith MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately 11620c0fd78eSBarry Smith after MatCreateShell() 11630c0fd78eSBarry Smith 11640c0fd78eSBarry Smith Logically Collective on Mat 11650c0fd78eSBarry Smith 11660c0fd78eSBarry Smith Input Parameter: 11670c0fd78eSBarry Smith . mat - the shell matrix 11680c0fd78eSBarry Smith 11690c0fd78eSBarry Smith Level: advanced 11700c0fd78eSBarry Smith 11710c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation() 11720c0fd78eSBarry Smith @*/ 11730c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A) 11740c0fd78eSBarry Smith { 11750c0fd78eSBarry Smith PetscErrorCode ierr; 1176976e8c5aSLisandro Dalcin Mat_Shell *shell; 11770c0fd78eSBarry Smith PetscBool flg; 11780c0fd78eSBarry Smith 11790c0fd78eSBarry Smith PetscFunctionBegin; 11800c0fd78eSBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 11810c0fd78eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr); 11820c0fd78eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 1183976e8c5aSLisandro Dalcin shell = (Mat_Shell*)A->data; 11840c0fd78eSBarry Smith shell->managescalingshifts = PETSC_FALSE; 1185976e8c5aSLisandro Dalcin A->ops->diagonalset = NULL; 1186976e8c5aSLisandro Dalcin A->ops->diagonalscale = NULL; 1187976e8c5aSLisandro Dalcin A->ops->scale = NULL; 1188976e8c5aSLisandro Dalcin A->ops->shift = NULL; 1189976e8c5aSLisandro Dalcin A->ops->axpy = NULL; 11900c0fd78eSBarry Smith PetscFunctionReturn(0); 11910c0fd78eSBarry Smith } 11920c0fd78eSBarry Smith 1193c16cb8f2SBarry Smith /*@C 1194f3b1f45cSBarry Smith MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function. 1195f3b1f45cSBarry Smith 1196f3b1f45cSBarry Smith Logically Collective on Mat 1197f3b1f45cSBarry Smith 1198f3b1f45cSBarry Smith Input Parameters: 1199f3b1f45cSBarry Smith + mat - the shell matrix 1200f3b1f45cSBarry Smith . f - the function 1201f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1202f3b1f45cSBarry Smith - ctx - an optional context for the function 1203f3b1f45cSBarry Smith 1204f3b1f45cSBarry Smith Output Parameter: 1205f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1206f3b1f45cSBarry Smith 1207f3b1f45cSBarry Smith Options Database: 1208f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1209f3b1f45cSBarry Smith 1210f3b1f45cSBarry Smith Level: advanced 1211f3b1f45cSBarry Smith 121295452b02SPatrick Sanan Fortran Notes: 121395452b02SPatrick Sanan Not supported from Fortran 1214f3b1f45cSBarry Smith 1215f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose() 1216f3b1f45cSBarry Smith @*/ 1217f3b1f45cSBarry Smith PetscErrorCode MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1218f3b1f45cSBarry Smith { 1219f3b1f45cSBarry Smith PetscErrorCode ierr; 1220f3b1f45cSBarry Smith PetscInt m,n; 1221f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1222f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 122374e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1224f3b1f45cSBarry Smith 1225f3b1f45cSBarry Smith PetscFunctionBegin; 1226f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1227f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr); 1228f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 1229f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr); 1230f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 1231f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 1232f3b1f45cSBarry Smith 12330bacdadaSStefano Zampini ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr); 12340bacdadaSStefano Zampini ierr = MatComputeOperator(mat,MATAIJ,&Dmat);CHKERRQ(ierr); 1235f3b1f45cSBarry Smith 1236f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 1237f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 1238f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 1239f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 1240f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1241f3b1f45cSBarry Smith flag = PETSC_FALSE; 1242f3b1f45cSBarry Smith if (v) { 1243fc7aafd1SBarry 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); 1244f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 1245f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 1246f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 1247f3b1f45cSBarry Smith } 1248f3b1f45cSBarry Smith } else if (v) { 1249fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr); 1250f3b1f45cSBarry Smith } 1251f3b1f45cSBarry Smith if (flg) *flg = flag; 1252f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 1253f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 1254f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 1255f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 1256f3b1f45cSBarry Smith PetscFunctionReturn(0); 1257f3b1f45cSBarry Smith } 1258f3b1f45cSBarry Smith 1259f3b1f45cSBarry Smith /*@C 1260f3b1f45cSBarry Smith MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function. 1261f3b1f45cSBarry Smith 1262f3b1f45cSBarry Smith Logically Collective on Mat 1263f3b1f45cSBarry Smith 1264f3b1f45cSBarry Smith Input Parameters: 1265f3b1f45cSBarry Smith + mat - the shell matrix 1266f3b1f45cSBarry Smith . f - the function 1267f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1268f3b1f45cSBarry Smith - ctx - an optional context for the function 1269f3b1f45cSBarry Smith 1270f3b1f45cSBarry Smith Output Parameter: 1271f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1272f3b1f45cSBarry Smith 1273f3b1f45cSBarry Smith Options Database: 1274f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1275f3b1f45cSBarry Smith 1276f3b1f45cSBarry Smith Level: advanced 1277f3b1f45cSBarry Smith 127895452b02SPatrick Sanan Fortran Notes: 127995452b02SPatrick Sanan Not supported from Fortran 1280f3b1f45cSBarry Smith 1281f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult() 1282f3b1f45cSBarry Smith @*/ 1283f3b1f45cSBarry Smith PetscErrorCode MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1284f3b1f45cSBarry Smith { 1285f3b1f45cSBarry Smith PetscErrorCode ierr; 1286f3b1f45cSBarry Smith Vec x,y,z; 1287f3b1f45cSBarry Smith PetscInt m,n,M,N; 1288f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1289f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 129074e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1291f3b1f45cSBarry Smith 1292f3b1f45cSBarry Smith PetscFunctionBegin; 1293f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1294f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr); 1295f3b1f45cSBarry Smith ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr); 1296f3b1f45cSBarry Smith ierr = VecDuplicate(y,&z);CHKERRQ(ierr); 1297f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 1298f3b1f45cSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 1299f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr); 1300f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 1301f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 13020bacdadaSStefano Zampini ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr); 1303f3b1f45cSBarry Smith ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr); 13040bacdadaSStefano Zampini ierr = MatComputeOperatorTranspose(mat,MATAIJ,&Dmat);CHKERRQ(ierr); 1305f3b1f45cSBarry Smith 1306f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 1307f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 1308f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 1309f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 1310f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1311f3b1f45cSBarry Smith flag = PETSC_FALSE; 1312f3b1f45cSBarry Smith if (v) { 1313fc7aafd1SBarry 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); 1314f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 1315f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 1316f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 1317f3b1f45cSBarry Smith } 1318f3b1f45cSBarry Smith } else if (v) { 1319fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr); 1320f3b1f45cSBarry Smith } 1321f3b1f45cSBarry Smith if (flg) *flg = flag; 1322f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 1323f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 1324f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 1325f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 1326f3b1f45cSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1327f3b1f45cSBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 1328f3b1f45cSBarry Smith ierr = VecDestroy(&z);CHKERRQ(ierr); 1329f3b1f45cSBarry Smith PetscFunctionReturn(0); 1330f3b1f45cSBarry Smith } 1331f3b1f45cSBarry Smith 1332f3b1f45cSBarry Smith /*@C 13330c0fd78eSBarry Smith MatShellSetOperation - Allows user to set a matrix operation for a shell matrix. 1334e51e0e81SBarry Smith 13353f9fe445SBarry Smith Logically Collective on Mat 1336fee21e36SBarry Smith 1337c7fcc2eaSBarry Smith Input Parameters: 1338c7fcc2eaSBarry Smith + mat - the shell matrix 1339c7fcc2eaSBarry Smith . op - the name of the operation 1340c7fcc2eaSBarry Smith - f - the function that provides the operation. 1341c7fcc2eaSBarry Smith 134215091d37SBarry Smith Level: advanced 134315091d37SBarry Smith 1344fae171e0SBarry Smith Usage: 1345e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 1346f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 1347c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 13480b627109SLois Curfman McInnes 1349a62d957aSLois Curfman McInnes Notes: 1350e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 13511c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 1352a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 13531c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 1354a62d957aSLois Curfman McInnes 1355a4d85fd7SBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 1356deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 1357deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 1358deebb3c3SLois Curfman McInnes routines, e.g., 1359a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1360a62d957aSLois Curfman McInnes 13614aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 13624aa34b0aSBarry Smith nonzero on failure. 13634aa34b0aSBarry Smith 1364a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 1365a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 1366a62d957aSLois Curfman McInnes set by MatCreateShell(). 1367a62d957aSLois Curfman McInnes 136895452b02SPatrick Sanan Fortran Notes: 136995452b02SPatrick Sanan For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 1370501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 1371501d9185SBarry Smith 13720c0fd78eSBarry Smith Use MatSetOperation() to set an operation for any matrix type 13730c0fd78eSBarry Smith 13740c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts() 1375e51e0e81SBarry Smith @*/ 13767087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 1377e51e0e81SBarry Smith { 1378ace3abfcSBarry Smith PetscBool flg; 1379976e8c5aSLisandro Dalcin Mat_Shell *shell; 1380976e8c5aSLisandro Dalcin PetscErrorCode ierr; 1381273d9f13SBarry Smith 13823a40ed3dSBarry Smith PetscFunctionBegin; 13830700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 13840c0fd78eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 13850c0fd78eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 1386976e8c5aSLisandro Dalcin shell = (Mat_Shell*)mat->data; 1387976e8c5aSLisandro Dalcin 13885edf6516SJed Brown switch (op) { 13895edf6516SJed Brown case MATOP_DESTROY: 139028f357e3SAlex Fikl shell->ops->destroy = (PetscErrorCode (*)(Mat))f; 13915edf6516SJed Brown break; 1392976e8c5aSLisandro Dalcin case MATOP_VIEW: 1393976e8c5aSLisandro Dalcin if (!mat->ops->viewnative) { 1394976e8c5aSLisandro Dalcin mat->ops->viewnative = mat->ops->view; 1395976e8c5aSLisandro Dalcin } 1396976e8c5aSLisandro Dalcin mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 1397976e8c5aSLisandro Dalcin break; 1398976e8c5aSLisandro Dalcin case MATOP_COPY: 1399976e8c5aSLisandro Dalcin shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f; 1400976e8c5aSLisandro Dalcin break; 14016464bf51SAlex Fikl case MATOP_DIAGONAL_SET: 14020c0fd78eSBarry Smith case MATOP_DIAGONAL_SCALE: 14030c0fd78eSBarry Smith case MATOP_SHIFT: 14040c0fd78eSBarry Smith case MATOP_SCALE: 14059f137db4SBarry Smith case MATOP_AXPY: 140645960306SStefano Zampini case MATOP_ZERO_ROWS: 140745960306SStefano Zampini case MATOP_ZERO_ROWS_COLUMNS: 14080c0fd78eSBarry Smith if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()"); 14090c0fd78eSBarry Smith (((void(**)(void))mat->ops)[op]) = f; 14106464bf51SAlex Fikl break; 14110c0fd78eSBarry Smith case MATOP_GET_DIAGONAL: 1412976e8c5aSLisandro Dalcin if (shell->managescalingshifts) { 1413976e8c5aSLisandro Dalcin shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1414976e8c5aSLisandro Dalcin mat->ops->getdiagonal = MatGetDiagonal_Shell; 1415976e8c5aSLisandro Dalcin } else { 1416976e8c5aSLisandro Dalcin shell->ops->getdiagonal = NULL; 1417976e8c5aSLisandro Dalcin mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 141840a5a6c4SBarry Smith } 14195edf6516SJed Brown break; 14205edf6516SJed Brown case MATOP_MULT: 14219f137db4SBarry Smith if (shell->managescalingshifts) { 14229f137db4SBarry Smith shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 14239f137db4SBarry Smith mat->ops->mult = MatMult_Shell; 1424976e8c5aSLisandro Dalcin } else { 1425976e8c5aSLisandro Dalcin shell->ops->mult = NULL; 1426976e8c5aSLisandro Dalcin mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1427976e8c5aSLisandro Dalcin } 14285edf6516SJed Brown break; 14295edf6516SJed Brown case MATOP_MULT_TRANSPOSE: 14309f137db4SBarry Smith if (shell->managescalingshifts) { 14319f137db4SBarry Smith shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 14325259c5a4SBarry Smith mat->ops->multtranspose = MatMultTranspose_Shell; 1433976e8c5aSLisandro Dalcin } else { 1434976e8c5aSLisandro Dalcin shell->ops->multtranspose = NULL; 1435976e8c5aSLisandro Dalcin mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1436976e8c5aSLisandro Dalcin } 14375edf6516SJed Brown break; 14385edf6516SJed Brown default: 14395edf6516SJed Brown (((void(**)(void))mat->ops)[op]) = f; 14405ab264f3SAlp Dener break; 1441a62d957aSLois Curfman McInnes } 14423a40ed3dSBarry Smith PetscFunctionReturn(0); 1443e51e0e81SBarry Smith } 1444f0479e8cSBarry Smith 1445d4bb536fSBarry Smith /*@C 1446d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 1447d4bb536fSBarry Smith 1448c7fcc2eaSBarry Smith Not Collective 1449c7fcc2eaSBarry Smith 1450d4bb536fSBarry Smith Input Parameters: 1451c7fcc2eaSBarry Smith + mat - the shell matrix 1452c7fcc2eaSBarry Smith - op - the name of the operation 1453d4bb536fSBarry Smith 1454d4bb536fSBarry Smith Output Parameter: 1455d4bb536fSBarry Smith . f - the function that provides the operation. 1456d4bb536fSBarry Smith 145715091d37SBarry Smith Level: advanced 145815091d37SBarry Smith 1459d4bb536fSBarry Smith Notes: 1460e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 1461d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 1462d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 1463d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 1464d4bb536fSBarry Smith 1465d4bb536fSBarry Smith All user-provided functions have the same calling 1466d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 1467d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 1468d4bb536fSBarry Smith routines, e.g., 1469d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1470d4bb536fSBarry Smith 1471d4bb536fSBarry Smith Within each user-defined routine, the user should call 1472d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 1473d4bb536fSBarry Smith set by MatCreateShell(). 1474d4bb536fSBarry Smith 1475ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 1476d4bb536fSBarry Smith @*/ 14777087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 1478d4bb536fSBarry Smith { 1479ace3abfcSBarry Smith PetscBool flg; 1480976e8c5aSLisandro Dalcin Mat_Shell *shell; 1481976e8c5aSLisandro Dalcin PetscErrorCode ierr; 1482273d9f13SBarry Smith 14833a40ed3dSBarry Smith PetscFunctionBegin; 14840700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1485976e8c5aSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 1486976e8c5aSLisandro Dalcin if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 1487976e8c5aSLisandro Dalcin shell = (Mat_Shell*)mat->data; 1488976e8c5aSLisandro Dalcin 148928f357e3SAlex Fikl switch (op) { 149028f357e3SAlex Fikl case MATOP_DESTROY: 149128f357e3SAlex Fikl *f = (void (*)(void))shell->ops->destroy; 149228f357e3SAlex Fikl break; 149328f357e3SAlex Fikl case MATOP_VIEW: 149428f357e3SAlex Fikl *f = (void (*)(void))mat->ops->view; 149528f357e3SAlex Fikl break; 1496976e8c5aSLisandro Dalcin case MATOP_COPY: 1497976e8c5aSLisandro Dalcin *f = (void (*)(void))shell->ops->copy; 1498976e8c5aSLisandro Dalcin break; 1499976e8c5aSLisandro Dalcin case MATOP_DIAGONAL_SET: 1500976e8c5aSLisandro Dalcin case MATOP_DIAGONAL_SCALE: 1501976e8c5aSLisandro Dalcin case MATOP_SHIFT: 1502976e8c5aSLisandro Dalcin case MATOP_SCALE: 1503976e8c5aSLisandro Dalcin case MATOP_AXPY: 150445960306SStefano Zampini case MATOP_ZERO_ROWS: 150545960306SStefano Zampini case MATOP_ZERO_ROWS_COLUMNS: 1506976e8c5aSLisandro Dalcin *f = (((void (**)(void))mat->ops)[op]); 1507976e8c5aSLisandro Dalcin break; 1508976e8c5aSLisandro Dalcin case MATOP_GET_DIAGONAL: 1509976e8c5aSLisandro Dalcin if (shell->ops->getdiagonal) 1510976e8c5aSLisandro Dalcin *f = (void (*)(void))shell->ops->getdiagonal; 1511976e8c5aSLisandro Dalcin else 1512976e8c5aSLisandro Dalcin *f = (((void (**)(void))mat->ops)[op]); 1513976e8c5aSLisandro Dalcin break; 1514976e8c5aSLisandro Dalcin case MATOP_MULT: 1515976e8c5aSLisandro Dalcin if (shell->ops->mult) 1516976e8c5aSLisandro Dalcin *f = (void (*)(void))shell->ops->mult; 1517976e8c5aSLisandro Dalcin else 1518976e8c5aSLisandro Dalcin *f = (((void (**)(void))mat->ops)[op]); 1519976e8c5aSLisandro Dalcin break; 1520976e8c5aSLisandro Dalcin case MATOP_MULT_TRANSPOSE: 1521976e8c5aSLisandro Dalcin if (shell->ops->multtranspose) 1522976e8c5aSLisandro Dalcin *f = (void (*)(void))shell->ops->multtranspose; 1523976e8c5aSLisandro Dalcin else 1524976e8c5aSLisandro Dalcin *f = (((void (**)(void))mat->ops)[op]); 1525976e8c5aSLisandro Dalcin break; 152628f357e3SAlex Fikl default: 1527c134de8dSSatish Balay *f = (((void (**)(void))mat->ops)[op]); 1528d4bb536fSBarry Smith } 15293a40ed3dSBarry Smith PetscFunctionReturn(0); 1530d4bb536fSBarry Smith } 1531