1be1d678aSKris Buschelman 2e51e0e81SBarry Smith /* 320563c6bSBarry Smith This provides a simple shell for Fortran (and C programmers) to 420563c6bSBarry Smith create a very simple matrix class for use with KSP without coding 5ed3cc1f0SBarry Smith much of anything. 6e51e0e81SBarry Smith */ 7e51e0e81SBarry Smith 8af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 945960306SStefano Zampini #include <petsc/private/vecimpl.h> 10e51e0e81SBarry Smith 1128f357e3SAlex Fikl struct _MatShellOps { 12976e8c5aSLisandro Dalcin /* 3 */ PetscErrorCode (*mult)(Mat, Vec, Vec); 13976e8c5aSLisandro Dalcin /* 5 */ PetscErrorCode (*multtranspose)(Mat, Vec, Vec); 14976e8c5aSLisandro Dalcin /* 17 */ PetscErrorCode (*getdiagonal)(Mat, Vec); 15976e8c5aSLisandro Dalcin /* 43 */ PetscErrorCode (*copy)(Mat, Mat, MatStructure); 16976e8c5aSLisandro Dalcin /* 60 */ PetscErrorCode (*destroy)(Mat); 1728f357e3SAlex Fikl }; 1828f357e3SAlex Fikl 19b77ba244SStefano Zampini struct _n_MatShellMatFunctionList { 20b77ba244SStefano Zampini PetscErrorCode (*symbolic)(Mat, Mat, Mat, void **); 21b77ba244SStefano Zampini PetscErrorCode (*numeric)(Mat, Mat, Mat, void *); 22b77ba244SStefano Zampini PetscErrorCode (*destroy)(void *); 23b77ba244SStefano Zampini MatProductType ptype; 24b77ba244SStefano Zampini char *composedname; /* string to identify routine with double dispatch */ 25b77ba244SStefano Zampini char *resultname; /* result matrix type */ 26b77ba244SStefano Zampini 27b77ba244SStefano Zampini struct _n_MatShellMatFunctionList *next; 28b77ba244SStefano Zampini }; 29b77ba244SStefano Zampini typedef struct _n_MatShellMatFunctionList *MatShellMatFunctionList; 30b77ba244SStefano Zampini 3128f357e3SAlex Fikl typedef struct { 3228f357e3SAlex Fikl struct _MatShellOps ops[1]; 332205254eSKarl Rupp 34b77ba244SStefano Zampini /* The user will manage the scaling and shifts for the MATSHELL, not the default */ 35b77ba244SStefano Zampini PetscBool managescalingshifts; 36b77ba244SStefano Zampini 37b77ba244SStefano Zampini /* support for MatScale, MatShift and MatMultAdd */ 38ef66eb69SBarry Smith PetscScalar vscale, vshift; 398fe8eb27SJed Brown Vec dshift; 408fe8eb27SJed Brown Vec left, right; 418fe8eb27SJed Brown Vec left_work, right_work; 425edf6516SJed Brown Vec left_add_work, right_add_work; 43b77ba244SStefano Zampini 44ef5c7bd2SStefano Zampini /* support for MatAXPY */ 459f137db4SBarry Smith Mat axpy; 469f137db4SBarry Smith PetscScalar axpy_vscale; 47b77ba244SStefano Zampini Vec axpy_left, axpy_right; 48ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 49b77ba244SStefano Zampini 5045960306SStefano Zampini /* support for ZeroRows/Columns operations */ 5145960306SStefano Zampini IS zrows; 5245960306SStefano Zampini IS zcols; 5345960306SStefano Zampini Vec zvals; 5445960306SStefano Zampini Vec zvals_w; 5545960306SStefano Zampini VecScatter zvals_sct_r; 5645960306SStefano Zampini VecScatter zvals_sct_c; 57b77ba244SStefano Zampini 58b77ba244SStefano Zampini /* MatMat operations */ 59b77ba244SStefano Zampini MatShellMatFunctionList matmat; 60b77ba244SStefano Zampini 61b77ba244SStefano Zampini /* user defined context */ 62800f99ffSJeremy L Thompson PetscContainer ctxcontainer; 6388cf3e7dSBarry Smith } Mat_Shell; 640c0fd78eSBarry Smith 6545960306SStefano Zampini /* 6645960306SStefano Zampini Store and scale values on zeroed rows 6745960306SStefano Zampini xx = [x_1, 0], 0 on zeroed columns 6845960306SStefano Zampini */ 699371c9d4SSatish Balay static PetscErrorCode MatShellPreZeroRight(Mat A, Vec x, Vec *xx) { 7045960306SStefano Zampini Mat_Shell *shell = (Mat_Shell *)A->data; 7145960306SStefano Zampini 7245960306SStefano Zampini PetscFunctionBegin; 7345960306SStefano Zampini *xx = x; 7445960306SStefano Zampini if (shell->zrows) { 759566063dSJacob Faibussowitsch PetscCall(VecSet(shell->zvals_w, 0.0)); 769566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_c, x, shell->zvals_w, INSERT_VALUES, SCATTER_FORWARD)); 779566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_c, x, shell->zvals_w, INSERT_VALUES, SCATTER_FORWARD)); 789566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->zvals_w, shell->zvals_w, shell->zvals)); 7945960306SStefano Zampini } 8045960306SStefano Zampini if (shell->zcols) { 8148a46eb9SPierre Jolivet if (!shell->right_work) PetscCall(MatCreateVecs(A, &shell->right_work, NULL)); 829566063dSJacob Faibussowitsch PetscCall(VecCopy(x, shell->right_work)); 839566063dSJacob Faibussowitsch PetscCall(VecISSet(shell->right_work, shell->zcols, 0.0)); 8445960306SStefano Zampini *xx = shell->right_work; 8545960306SStefano Zampini } 8645960306SStefano Zampini PetscFunctionReturn(0); 8745960306SStefano Zampini } 8845960306SStefano Zampini 8945960306SStefano Zampini /* Insert properly diagonally scaled values stored in MatShellPreZeroRight */ 909371c9d4SSatish Balay static PetscErrorCode MatShellPostZeroLeft(Mat A, Vec x) { 9145960306SStefano Zampini Mat_Shell *shell = (Mat_Shell *)A->data; 9245960306SStefano Zampini 9345960306SStefano Zampini PetscFunctionBegin; 9445960306SStefano Zampini if (shell->zrows) { 959566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_r, shell->zvals_w, x, INSERT_VALUES, SCATTER_REVERSE)); 969566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_r, shell->zvals_w, x, INSERT_VALUES, SCATTER_REVERSE)); 9745960306SStefano Zampini } 9845960306SStefano Zampini PetscFunctionReturn(0); 9945960306SStefano Zampini } 10045960306SStefano Zampini 10145960306SStefano Zampini /* 10245960306SStefano Zampini Store and scale values on zeroed rows 10345960306SStefano Zampini xx = [x_1, 0], 0 on zeroed rows 10445960306SStefano Zampini */ 1059371c9d4SSatish Balay static PetscErrorCode MatShellPreZeroLeft(Mat A, Vec x, Vec *xx) { 10645960306SStefano Zampini Mat_Shell *shell = (Mat_Shell *)A->data; 10745960306SStefano Zampini 10845960306SStefano Zampini PetscFunctionBegin; 10945960306SStefano Zampini *xx = NULL; 11045960306SStefano Zampini if (!shell->zrows) { 11145960306SStefano Zampini *xx = x; 11245960306SStefano Zampini } else { 11348a46eb9SPierre Jolivet if (!shell->left_work) PetscCall(MatCreateVecs(A, NULL, &shell->left_work)); 1149566063dSJacob Faibussowitsch PetscCall(VecCopy(x, shell->left_work)); 1159566063dSJacob Faibussowitsch PetscCall(VecSet(shell->zvals_w, 0.0)); 1169566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_r, shell->zvals_w, shell->left_work, INSERT_VALUES, SCATTER_REVERSE)); 1179566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_r, shell->zvals_w, shell->left_work, INSERT_VALUES, SCATTER_REVERSE)); 1189566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_r, x, shell->zvals_w, INSERT_VALUES, SCATTER_FORWARD)); 1199566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_r, x, shell->zvals_w, INSERT_VALUES, SCATTER_FORWARD)); 1209566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->zvals_w, shell->zvals_w, shell->zvals)); 12145960306SStefano Zampini *xx = shell->left_work; 12245960306SStefano Zampini } 12345960306SStefano Zampini PetscFunctionReturn(0); 12445960306SStefano Zampini } 12545960306SStefano Zampini 12645960306SStefano Zampini /* Zero zero-columns contributions, sum contributions from properly scaled values stored in MatShellPreZeroLeft */ 1279371c9d4SSatish Balay static PetscErrorCode MatShellPostZeroRight(Mat A, Vec x) { 12845960306SStefano Zampini Mat_Shell *shell = (Mat_Shell *)A->data; 12945960306SStefano Zampini 13045960306SStefano Zampini PetscFunctionBegin; 1311baa6e33SBarry Smith if (shell->zcols) PetscCall(VecISSet(x, shell->zcols, 0.0)); 13245960306SStefano Zampini if (shell->zrows) { 1339566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_c, shell->zvals_w, x, ADD_VALUES, SCATTER_REVERSE)); 1349566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_c, shell->zvals_w, x, ADD_VALUES, SCATTER_REVERSE)); 13545960306SStefano Zampini } 13645960306SStefano Zampini PetscFunctionReturn(0); 13745960306SStefano Zampini } 13845960306SStefano Zampini 1398fe8eb27SJed Brown /* 1400c0fd78eSBarry Smith xx = diag(left)*x 1418fe8eb27SJed Brown */ 1429371c9d4SSatish Balay static PetscErrorCode MatShellPreScaleLeft(Mat A, Vec x, Vec *xx) { 1438fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell *)A->data; 1448fe8eb27SJed Brown 1458fe8eb27SJed Brown PetscFunctionBegin; 1460298fd71SBarry Smith *xx = NULL; 1478fe8eb27SJed Brown if (!shell->left) { 1488fe8eb27SJed Brown *xx = x; 1498fe8eb27SJed Brown } else { 1509566063dSJacob Faibussowitsch if (!shell->left_work) PetscCall(VecDuplicate(shell->left, &shell->left_work)); 1519566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->left_work, x, shell->left)); 1528fe8eb27SJed Brown *xx = shell->left_work; 1538fe8eb27SJed Brown } 1548fe8eb27SJed Brown PetscFunctionReturn(0); 1558fe8eb27SJed Brown } 1568fe8eb27SJed Brown 1570c0fd78eSBarry Smith /* 1580c0fd78eSBarry Smith xx = diag(right)*x 1590c0fd78eSBarry Smith */ 1609371c9d4SSatish Balay static PetscErrorCode MatShellPreScaleRight(Mat A, Vec x, Vec *xx) { 1618fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell *)A->data; 1628fe8eb27SJed Brown 1638fe8eb27SJed Brown PetscFunctionBegin; 1640298fd71SBarry Smith *xx = NULL; 1658fe8eb27SJed Brown if (!shell->right) { 1668fe8eb27SJed Brown *xx = x; 1678fe8eb27SJed Brown } else { 1689566063dSJacob Faibussowitsch if (!shell->right_work) PetscCall(VecDuplicate(shell->right, &shell->right_work)); 1699566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->right_work, x, shell->right)); 1708fe8eb27SJed Brown *xx = shell->right_work; 1718fe8eb27SJed Brown } 1728fe8eb27SJed Brown PetscFunctionReturn(0); 1738fe8eb27SJed Brown } 1748fe8eb27SJed Brown 1750c0fd78eSBarry Smith /* 1760c0fd78eSBarry Smith x = diag(left)*x 1770c0fd78eSBarry Smith */ 1789371c9d4SSatish Balay static PetscErrorCode MatShellPostScaleLeft(Mat A, Vec x) { 1798fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell *)A->data; 1808fe8eb27SJed Brown 1818fe8eb27SJed Brown PetscFunctionBegin; 1829566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseMult(x, x, shell->left)); 1838fe8eb27SJed Brown PetscFunctionReturn(0); 1848fe8eb27SJed Brown } 1858fe8eb27SJed Brown 1860c0fd78eSBarry Smith /* 1870c0fd78eSBarry Smith x = diag(right)*x 1880c0fd78eSBarry Smith */ 1899371c9d4SSatish Balay static PetscErrorCode MatShellPostScaleRight(Mat A, Vec x) { 1908fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell *)A->data; 1918fe8eb27SJed Brown 1928fe8eb27SJed Brown PetscFunctionBegin; 1939566063dSJacob Faibussowitsch if (shell->right) PetscCall(VecPointwiseMult(x, x, shell->right)); 1948fe8eb27SJed Brown PetscFunctionReturn(0); 1958fe8eb27SJed Brown } 1968fe8eb27SJed Brown 1970c0fd78eSBarry Smith /* 1980c0fd78eSBarry Smith Y = vscale*Y + diag(dshift)*X + vshift*X 1990c0fd78eSBarry Smith 2000c0fd78eSBarry Smith On input Y already contains A*x 2010c0fd78eSBarry Smith */ 2029371c9d4SSatish Balay static PetscErrorCode MatShellShiftAndScale(Mat A, Vec X, Vec Y) { 2038fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell *)A->data; 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; 2109566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &m)); 2119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(shell->dshift, &d)); 2129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2139566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, &y)); 2148fe8eb27SJed Brown for (i = 0; i < m; i++) y[i] = shell->vscale * y[i] + d[i] * x[i]; 2159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(shell->dshift, &d)); 2169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, &y)); 2180c0fd78eSBarry Smith } else { 2199566063dSJacob Faibussowitsch PetscCall(VecScale(Y, shell->vscale)); 2208fe8eb27SJed Brown } 2219566063dSJacob Faibussowitsch if (shell->vshift != 0.0) PetscCall(VecAXPY(Y, shell->vshift, X)); /* if test is for non-square matrices */ 2228fe8eb27SJed Brown PetscFunctionReturn(0); 2238fe8eb27SJed Brown } 224e51e0e81SBarry Smith 2259371c9d4SSatish Balay static PetscErrorCode MatShellGetContext_Shell(Mat mat, void *ctx) { 226800f99ffSJeremy L Thompson Mat_Shell *shell = (Mat_Shell *)mat->data; 227800f99ffSJeremy L Thompson 228789d8953SBarry Smith PetscFunctionBegin; 229800f99ffSJeremy L Thompson if (shell->ctxcontainer) PetscCall(PetscContainerGetPointer(shell->ctxcontainer, (void **)ctx)); 230800f99ffSJeremy L Thompson else *(void **)ctx = NULL; 231789d8953SBarry Smith PetscFunctionReturn(0); 232789d8953SBarry Smith } 233789d8953SBarry Smith 2349d225801SJed Brown /*@ 23511a5261eSBarry Smith MatShellGetContext - Returns the user-provided context associated with a `MATSHELL` shell matrix. 236b4fd4287SBarry Smith 237c7fcc2eaSBarry Smith Not Collective 238c7fcc2eaSBarry Smith 239b4fd4287SBarry Smith Input Parameter: 24011a5261eSBarry Smith . mat - the matrix, should have been created with `MatCreateShell()` 241b4fd4287SBarry Smith 242b4fd4287SBarry Smith Output Parameter: 243b4fd4287SBarry Smith . ctx - the user provided context 244b4fd4287SBarry Smith 24515091d37SBarry Smith Level: advanced 24615091d37SBarry Smith 24711a5261eSBarry Smith Fortran Note: 24895452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 249daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 250a62d957aSLois Curfman McInnes 25111a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateShell()`, `MatShellSetOperation()`, `MatShellSetContext()` 2520bc0a719SBarry Smith @*/ 2539371c9d4SSatish Balay PetscErrorCode MatShellGetContext(Mat mat, void *ctx) { 2543a40ed3dSBarry Smith PetscFunctionBegin; 2550700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2564482741eSBarry Smith PetscValidPointer(ctx, 2); 257cac4c232SBarry Smith PetscUseMethod(mat, "MatShellGetContext_C", (Mat, void *), (mat, ctx)); 2583a40ed3dSBarry Smith PetscFunctionReturn(0); 259b4fd4287SBarry Smith } 260b4fd4287SBarry Smith 2619371c9d4SSatish Balay static PetscErrorCode MatZeroRowsColumns_Local_Shell(Mat mat, PetscInt nr, PetscInt rows[], PetscInt nc, PetscInt cols[], PetscScalar diag, PetscBool rc) { 26245960306SStefano Zampini Mat_Shell *shell = (Mat_Shell *)mat->data; 26345960306SStefano Zampini Vec x = NULL, b = NULL; 26445960306SStefano Zampini IS is1, is2; 26545960306SStefano Zampini const PetscInt *ridxs; 26645960306SStefano Zampini PetscInt *idxs, *gidxs; 26745960306SStefano Zampini PetscInt cum, rst, cst, i; 26845960306SStefano Zampini 26945960306SStefano Zampini PetscFunctionBegin; 27048a46eb9SPierre Jolivet if (!shell->zvals) PetscCall(MatCreateVecs(mat, NULL, &shell->zvals)); 27148a46eb9SPierre Jolivet if (!shell->zvals_w) PetscCall(VecDuplicate(shell->zvals, &shell->zvals_w)); 2729566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat, &rst, NULL)); 2739566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(mat, &cst, NULL)); 27445960306SStefano Zampini 27545960306SStefano Zampini /* Expand/create index set of zeroed rows */ 2769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nr, &idxs)); 27745960306SStefano Zampini for (i = 0; i < nr; i++) idxs[i] = rows[i] + rst; 2789566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nr, idxs, PETSC_OWN_POINTER, &is1)); 2799566063dSJacob Faibussowitsch PetscCall(ISSort(is1)); 2809566063dSJacob Faibussowitsch PetscCall(VecISSet(shell->zvals, is1, diag)); 28145960306SStefano Zampini if (shell->zrows) { 2829566063dSJacob Faibussowitsch PetscCall(ISSum(shell->zrows, is1, &is2)); 2839566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zrows)); 2849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 28545960306SStefano Zampini shell->zrows = is2; 28645960306SStefano Zampini } else shell->zrows = is1; 28745960306SStefano Zampini 28845960306SStefano Zampini /* Create scatters for diagonal values communications */ 2899566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_c)); 2909566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_r)); 29145960306SStefano Zampini 29245960306SStefano Zampini /* row scatter: from/to left vector */ 2939566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat, &x, &b)); 2949566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(b, shell->zrows, shell->zvals_w, shell->zrows, &shell->zvals_sct_r)); 29545960306SStefano Zampini 29645960306SStefano Zampini /* col scatter: from right vector to left vector */ 2979566063dSJacob Faibussowitsch PetscCall(ISGetIndices(shell->zrows, &ridxs)); 2989566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(shell->zrows, &nr)); 2999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nr, &gidxs)); 30045960306SStefano Zampini for (i = 0, cum = 0; i < nr; i++) { 30145960306SStefano Zampini if (ridxs[i] >= mat->cmap->N) continue; 30245960306SStefano Zampini gidxs[cum] = ridxs[i]; 30345960306SStefano Zampini cum++; 30445960306SStefano Zampini } 3059566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)mat), cum, gidxs, PETSC_OWN_POINTER, &is1)); 3069566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, is1, shell->zvals_w, is1, &shell->zvals_sct_c)); 3079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 3089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 31045960306SStefano Zampini 31145960306SStefano Zampini /* Expand/create index set of zeroed columns */ 31245960306SStefano Zampini if (rc) { 3139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nc, &idxs)); 31445960306SStefano Zampini for (i = 0; i < nc; i++) idxs[i] = cols[i] + cst; 3159566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nc, idxs, PETSC_OWN_POINTER, &is1)); 3169566063dSJacob Faibussowitsch PetscCall(ISSort(is1)); 31745960306SStefano Zampini if (shell->zcols) { 3189566063dSJacob Faibussowitsch PetscCall(ISSum(shell->zcols, is1, &is2)); 3199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zcols)); 3209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 32145960306SStefano Zampini shell->zcols = is2; 32245960306SStefano Zampini } else shell->zcols = is1; 32345960306SStefano Zampini } 32445960306SStefano Zampini PetscFunctionReturn(0); 32545960306SStefano Zampini } 32645960306SStefano Zampini 3279371c9d4SSatish Balay static PetscErrorCode MatZeroRows_Shell(Mat mat, PetscInt n, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) { 328ef5c7bd2SStefano Zampini Mat_Shell *shell = (Mat_Shell *)mat->data; 32945960306SStefano Zampini PetscInt nr, *lrows; 33045960306SStefano Zampini 33145960306SStefano Zampini PetscFunctionBegin; 33245960306SStefano Zampini if (x && b) { 33345960306SStefano Zampini Vec xt; 33445960306SStefano Zampini PetscScalar *vals; 33545960306SStefano Zampini PetscInt *gcols, i, st, nl, nc; 33645960306SStefano Zampini 3379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &gcols)); 3389371c9d4SSatish Balay for (i = 0, nc = 0; i < n; i++) 3399371c9d4SSatish Balay if (rows[i] < mat->cmap->N) gcols[nc++] = rows[i]; 34045960306SStefano Zampini 3419566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat, &xt, NULL)); 3429566063dSJacob Faibussowitsch PetscCall(VecCopy(x, xt)); 3439566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nc, &vals)); 3449566063dSJacob Faibussowitsch PetscCall(VecSetValues(xt, nc, gcols, vals, INSERT_VALUES)); /* xt = [x1, 0] */ 3459566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 3469566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(xt)); 3479566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(xt)); 3489566063dSJacob Faibussowitsch PetscCall(VecAYPX(xt, -1.0, x)); /* xt = [0, x2] */ 34945960306SStefano Zampini 3509566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(xt, &st, NULL)); 3519566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xt, &nl)); 3529566063dSJacob Faibussowitsch PetscCall(VecGetArray(xt, &vals)); 35345960306SStefano Zampini for (i = 0; i < nl; i++) { 35445960306SStefano Zampini PetscInt g = i + st; 35545960306SStefano Zampini if (g > mat->rmap->N) continue; 35645960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 3579566063dSJacob Faibussowitsch PetscCall(VecSetValue(b, g, diag * vals[i], INSERT_VALUES)); 35845960306SStefano Zampini } 3599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xt, &vals)); 3609566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b)); 3619566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b)); /* b = [b1, x2 * diag] */ 3629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xt)); 3639566063dSJacob Faibussowitsch PetscCall(PetscFree(gcols)); 36445960306SStefano Zampini } 3659566063dSJacob Faibussowitsch PetscCall(PetscLayoutMapLocal(mat->rmap, n, rows, &nr, &lrows, NULL)); 3669566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns_Local_Shell(mat, nr, lrows, 0, NULL, diag, PETSC_FALSE)); 3671baa6e33SBarry Smith if (shell->axpy) PetscCall(MatZeroRows(shell->axpy, n, rows, 0.0, NULL, NULL)); 3689566063dSJacob Faibussowitsch PetscCall(PetscFree(lrows)); 36945960306SStefano Zampini PetscFunctionReturn(0); 37045960306SStefano Zampini } 37145960306SStefano Zampini 3729371c9d4SSatish Balay static PetscErrorCode MatZeroRowsColumns_Shell(Mat mat, PetscInt n, const PetscInt rowscols[], PetscScalar diag, Vec x, Vec b) { 373ef5c7bd2SStefano Zampini Mat_Shell *shell = (Mat_Shell *)mat->data; 37445960306SStefano Zampini PetscInt *lrows, *lcols; 37545960306SStefano Zampini PetscInt nr, nc; 37645960306SStefano Zampini PetscBool congruent; 37745960306SStefano Zampini 37845960306SStefano Zampini PetscFunctionBegin; 37945960306SStefano Zampini if (x && b) { 38045960306SStefano Zampini Vec xt, bt; 38145960306SStefano Zampini PetscScalar *vals; 38245960306SStefano Zampini PetscInt *grows, *gcols, i, st, nl; 38345960306SStefano Zampini 3849566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &grows, n, &gcols)); 3859371c9d4SSatish Balay for (i = 0, nr = 0; i < n; i++) 3869371c9d4SSatish Balay if (rowscols[i] < mat->rmap->N) grows[nr++] = rowscols[i]; 3879371c9d4SSatish Balay for (i = 0, nc = 0; i < n; i++) 3889371c9d4SSatish Balay if (rowscols[i] < mat->cmap->N) gcols[nc++] = rowscols[i]; 3899566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(n, &vals)); 39045960306SStefano Zampini 3919566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat, &xt, &bt)); 3929566063dSJacob Faibussowitsch PetscCall(VecCopy(x, xt)); 3939566063dSJacob Faibussowitsch PetscCall(VecSetValues(xt, nc, gcols, vals, INSERT_VALUES)); /* xt = [x1, 0] */ 3949566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(xt)); 3959566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(xt)); 3969566063dSJacob Faibussowitsch PetscCall(VecAXPY(xt, -1.0, x)); /* xt = [0, -x2] */ 3979566063dSJacob Faibussowitsch PetscCall(MatMult(mat, xt, bt)); /* bt = [-A12*x2,-A22*x2] */ 3989566063dSJacob Faibussowitsch PetscCall(VecSetValues(bt, nr, grows, vals, INSERT_VALUES)); /* bt = [-A12*x2,0] */ 3999566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bt)); 4009566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bt)); 4019566063dSJacob Faibussowitsch PetscCall(VecAXPY(b, 1.0, bt)); /* b = [b1 - A12*x2, b2] */ 4029566063dSJacob Faibussowitsch PetscCall(VecSetValues(bt, nr, grows, vals, INSERT_VALUES)); /* b = [b1 - A12*x2, 0] */ 4039566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bt)); 4049566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bt)); 4059566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 40645960306SStefano Zampini 4079566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(xt, &st, NULL)); 4089566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xt, &nl)); 4099566063dSJacob Faibussowitsch PetscCall(VecGetArray(xt, &vals)); 41045960306SStefano Zampini for (i = 0; i < nl; i++) { 41145960306SStefano Zampini PetscInt g = i + st; 41245960306SStefano Zampini if (g > mat->rmap->N) continue; 41345960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 4149566063dSJacob Faibussowitsch PetscCall(VecSetValue(b, g, -diag * vals[i], INSERT_VALUES)); 41545960306SStefano Zampini } 4169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xt, &vals)); 4179566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b)); 4189566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b)); /* b = [b1 - A12*x2, x2 * diag] */ 4199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xt)); 4209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bt)); 4219566063dSJacob Faibussowitsch PetscCall(PetscFree2(grows, gcols)); 42245960306SStefano Zampini } 4239566063dSJacob Faibussowitsch PetscCall(PetscLayoutMapLocal(mat->rmap, n, rowscols, &nr, &lrows, NULL)); 4249566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(mat, &congruent)); 42545960306SStefano Zampini if (congruent) { 42645960306SStefano Zampini nc = nr; 42745960306SStefano Zampini lcols = lrows; 42845960306SStefano Zampini } else { /* MatZeroRowsColumns implicitly assumes the rowscols indices are for a square matrix, here we handle a more general case */ 42945960306SStefano Zampini PetscInt i, nt, *t; 43045960306SStefano Zampini 4319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &t)); 4329371c9d4SSatish Balay for (i = 0, nt = 0; i < n; i++) 4339371c9d4SSatish Balay if (rowscols[i] < mat->cmap->N) t[nt++] = rowscols[i]; 4349566063dSJacob Faibussowitsch PetscCall(PetscLayoutMapLocal(mat->cmap, nt, t, &nc, &lcols, NULL)); 4359566063dSJacob Faibussowitsch PetscCall(PetscFree(t)); 43645960306SStefano Zampini } 4379566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns_Local_Shell(mat, nr, lrows, nc, lcols, diag, PETSC_TRUE)); 43848a46eb9SPierre Jolivet if (!congruent) PetscCall(PetscFree(lcols)); 4399566063dSJacob Faibussowitsch PetscCall(PetscFree(lrows)); 4401baa6e33SBarry Smith if (shell->axpy) PetscCall(MatZeroRowsColumns(shell->axpy, n, rowscols, 0.0, NULL, NULL)); 44145960306SStefano Zampini PetscFunctionReturn(0); 44245960306SStefano Zampini } 44345960306SStefano Zampini 4449371c9d4SSatish Balay static PetscErrorCode MatDestroy_Shell(Mat mat) { 445bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell *)mat->data; 446b77ba244SStefano Zampini MatShellMatFunctionList matmat; 447ed3cc1f0SBarry Smith 4483a40ed3dSBarry Smith PetscFunctionBegin; 4491baa6e33SBarry Smith if (shell->ops->destroy) PetscCall((*shell->ops->destroy)(mat)); 4509566063dSJacob Faibussowitsch PetscCall(PetscMemzero(shell->ops, sizeof(struct _MatShellOps))); 4519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->left)); 4529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->right)); 4539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->dshift)); 4549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->left_work)); 4559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->right_work)); 4569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->left_add_work)); 4579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->right_add_work)); 4589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->axpy_left)); 4599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->axpy_right)); 4609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&shell->axpy)); 4619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->zvals_w)); 4629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->zvals)); 4639566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_c)); 4649566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_r)); 4659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zrows)); 4669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zcols)); 467b77ba244SStefano Zampini 468b77ba244SStefano Zampini matmat = shell->matmat; 469b77ba244SStefano Zampini while (matmat) { 470b77ba244SStefano Zampini MatShellMatFunctionList next = matmat->next; 471b77ba244SStefano Zampini 4729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, matmat->composedname, NULL)); 4739566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat->composedname)); 4749566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat->resultname)); 4759566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat)); 476b77ba244SStefano Zampini matmat = next; 477b77ba244SStefano Zampini } 478800f99ffSJeremy L Thompson PetscCall(MatShellSetContext(mat, NULL)); 4799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellGetContext_C", NULL)); 4809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetContext_C", NULL)); 481800f99ffSJeremy L Thompson PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetContextDestroy_C", NULL)); 4829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetVecType_C", NULL)); 4839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetManageScalingShifts_C", NULL)); 4849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetOperation_C", NULL)); 4859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellGetOperation_C", NULL)); 4869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetMatProductOperation_C", NULL)); 4879566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 488b77ba244SStefano Zampini PetscFunctionReturn(0); 489b77ba244SStefano Zampini } 490b77ba244SStefano Zampini 491b77ba244SStefano Zampini typedef struct { 492b77ba244SStefano Zampini PetscErrorCode (*numeric)(Mat, Mat, Mat, void *); 493b77ba244SStefano Zampini PetscErrorCode (*destroy)(void *); 494b77ba244SStefano Zampini void *userdata; 495b77ba244SStefano Zampini Mat B; 496b77ba244SStefano Zampini Mat Bt; 497b77ba244SStefano Zampini Mat axpy; 498b77ba244SStefano Zampini } MatMatDataShell; 499b77ba244SStefano Zampini 5009371c9d4SSatish Balay static PetscErrorCode DestroyMatMatDataShell(void *data) { 501b77ba244SStefano Zampini MatMatDataShell *mmdata = (MatMatDataShell *)data; 502b77ba244SStefano Zampini 503b77ba244SStefano Zampini PetscFunctionBegin; 5041baa6e33SBarry Smith if (mmdata->destroy) PetscCall((*mmdata->destroy)(mmdata->userdata)); 5059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mmdata->B)); 5069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mmdata->Bt)); 5079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mmdata->axpy)); 5089566063dSJacob Faibussowitsch PetscCall(PetscFree(mmdata)); 509b77ba244SStefano Zampini PetscFunctionReturn(0); 510b77ba244SStefano Zampini } 511b77ba244SStefano Zampini 5129371c9d4SSatish Balay static PetscErrorCode MatProductNumeric_Shell_X(Mat D) { 513b77ba244SStefano Zampini Mat_Product *product; 514b77ba244SStefano Zampini Mat A, B; 515b77ba244SStefano Zampini MatMatDataShell *mdata; 516b77ba244SStefano Zampini PetscScalar zero = 0.0; 517b77ba244SStefano Zampini 518b77ba244SStefano Zampini PetscFunctionBegin; 519b77ba244SStefano Zampini MatCheckProduct(D, 1); 520b77ba244SStefano Zampini product = D->product; 52128b400f6SJacob Faibussowitsch PetscCheck(product->data, PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "Product data empty"); 522b77ba244SStefano Zampini A = product->A; 523b77ba244SStefano Zampini B = product->B; 524b77ba244SStefano Zampini mdata = (MatMatDataShell *)product->data; 525b77ba244SStefano Zampini if (mdata->numeric) { 526b77ba244SStefano Zampini Mat_Shell *shell = (Mat_Shell *)A->data; 527b77ba244SStefano Zampini PetscErrorCode (*stashsym)(Mat) = D->ops->productsymbolic; 528b77ba244SStefano Zampini PetscErrorCode (*stashnum)(Mat) = D->ops->productnumeric; 529b77ba244SStefano Zampini PetscBool useBmdata = PETSC_FALSE, newB = PETSC_TRUE; 530b77ba244SStefano Zampini 531b77ba244SStefano Zampini if (shell->managescalingshifts) { 53208401ef6SPierre Jolivet PetscCheck(!shell->zcols && !shell->zrows, PetscObjectComm((PetscObject)D), PETSC_ERR_SUP, "MatProduct not supported with zeroed rows/columns"); 533b77ba244SStefano Zampini if (shell->right || shell->left) { 534b77ba244SStefano Zampini useBmdata = PETSC_TRUE; 535b77ba244SStefano Zampini if (!mdata->B) { 5369566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_SHARE_NONZERO_PATTERN, &mdata->B)); 537b77ba244SStefano Zampini } else { 538b77ba244SStefano Zampini newB = PETSC_FALSE; 539b77ba244SStefano Zampini } 5409566063dSJacob Faibussowitsch PetscCall(MatCopy(B, mdata->B, SAME_NONZERO_PATTERN)); 541b77ba244SStefano Zampini } 542b77ba244SStefano Zampini switch (product->type) { 543b77ba244SStefano Zampini case MATPRODUCT_AB: /* s L A R B + v L R B + L D R B */ 5441baa6e33SBarry Smith if (shell->right) PetscCall(MatDiagonalScale(mdata->B, shell->right, NULL)); 545b77ba244SStefano Zampini break; 546b77ba244SStefano Zampini case MATPRODUCT_AtB: /* s R A^t L B + v R L B + R D L B */ 5471baa6e33SBarry Smith if (shell->left) PetscCall(MatDiagonalScale(mdata->B, shell->left, NULL)); 548b77ba244SStefano Zampini break; 549b77ba244SStefano Zampini case MATPRODUCT_ABt: /* s L A R B^t + v L R B^t + L D R B^t */ 5501baa6e33SBarry Smith if (shell->right) PetscCall(MatDiagonalScale(mdata->B, NULL, shell->right)); 551b77ba244SStefano Zampini break; 552b77ba244SStefano Zampini case MATPRODUCT_RARt: /* s B L A R B^t + v B L R B^t + B L D R B^t */ 553b77ba244SStefano Zampini if (shell->right && shell->left) { 554b77ba244SStefano Zampini PetscBool flg; 555b77ba244SStefano Zampini 5569566063dSJacob Faibussowitsch PetscCall(VecEqual(shell->right, shell->left, &flg)); 5579371c9d4SSatish Balay PetscCheck(flg, PetscObjectComm((PetscObject)D), PETSC_ERR_SUP, "MatProductSymbolic type %s not supported for %s and %s matrices because left scaling != from right scaling", MatProductTypes[product->type], ((PetscObject)A)->type_name, 5589371c9d4SSatish Balay ((PetscObject)B)->type_name); 559b77ba244SStefano Zampini } 5601baa6e33SBarry Smith if (shell->right) PetscCall(MatDiagonalScale(mdata->B, NULL, shell->right)); 561b77ba244SStefano Zampini break; 562b77ba244SStefano Zampini case MATPRODUCT_PtAP: /* s B^t L A R B + v B^t L R B + B^t L D R B */ 563b77ba244SStefano Zampini if (shell->right && shell->left) { 564b77ba244SStefano Zampini PetscBool flg; 565b77ba244SStefano Zampini 5669566063dSJacob Faibussowitsch PetscCall(VecEqual(shell->right, shell->left, &flg)); 5679371c9d4SSatish Balay PetscCheck(flg, PetscObjectComm((PetscObject)D), PETSC_ERR_SUP, "MatProductSymbolic type %s not supported for %s and %s matrices because left scaling != from right scaling", MatProductTypes[product->type], ((PetscObject)A)->type_name, 5689371c9d4SSatish Balay ((PetscObject)B)->type_name); 569b77ba244SStefano Zampini } 5701baa6e33SBarry Smith if (shell->right) PetscCall(MatDiagonalScale(mdata->B, shell->right, NULL)); 571b77ba244SStefano Zampini break; 57298921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)D), PETSC_ERR_SUP, "MatProductSymbolic type %s not supported for %s and %s matrices", MatProductTypes[product->type], ((PetscObject)A)->type_name, ((PetscObject)B)->type_name); 573b77ba244SStefano Zampini } 574b77ba244SStefano Zampini } 575b77ba244SStefano Zampini /* allow the user to call MatMat operations on D */ 576b77ba244SStefano Zampini D->product = NULL; 577b77ba244SStefano Zampini D->ops->productsymbolic = NULL; 578b77ba244SStefano Zampini D->ops->productnumeric = NULL; 579b77ba244SStefano Zampini 5809566063dSJacob Faibussowitsch PetscCall((*mdata->numeric)(A, useBmdata ? mdata->B : B, D, mdata->userdata)); 581b77ba244SStefano Zampini 582b77ba244SStefano Zampini /* clear any leftover user data and restore D pointers */ 5839566063dSJacob Faibussowitsch PetscCall(MatProductClear(D)); 584b77ba244SStefano Zampini D->ops->productsymbolic = stashsym; 585b77ba244SStefano Zampini D->ops->productnumeric = stashnum; 586b77ba244SStefano Zampini D->product = product; 587b77ba244SStefano Zampini 588b77ba244SStefano Zampini if (shell->managescalingshifts) { 5899566063dSJacob Faibussowitsch PetscCall(MatScale(D, shell->vscale)); 590b77ba244SStefano Zampini switch (product->type) { 591b77ba244SStefano Zampini case MATPRODUCT_AB: /* s L A R B + v L R B + L D R B */ 592b77ba244SStefano Zampini case MATPRODUCT_ABt: /* s L A R B^t + v L R B^t + L D R B^t */ 593b77ba244SStefano Zampini if (shell->left) { 5949566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(D, shell->left, NULL)); 595b77ba244SStefano Zampini if (shell->dshift || shell->vshift != zero) { 5969566063dSJacob Faibussowitsch if (!shell->left_work) PetscCall(MatCreateVecs(A, NULL, &shell->left_work)); 597b77ba244SStefano Zampini if (shell->dshift) { 5989566063dSJacob Faibussowitsch PetscCall(VecCopy(shell->dshift, shell->left_work)); 5999566063dSJacob Faibussowitsch PetscCall(VecShift(shell->left_work, shell->vshift)); 6009566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->left_work, shell->left_work, shell->left)); 601b77ba244SStefano Zampini } else { 6029566063dSJacob Faibussowitsch PetscCall(VecSet(shell->left_work, shell->vshift)); 603b77ba244SStefano Zampini } 604b77ba244SStefano Zampini if (product->type == MATPRODUCT_ABt) { 605b77ba244SStefano Zampini MatReuse reuse = mdata->Bt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX; 606b77ba244SStefano Zampini MatStructure str = mdata->Bt ? SUBSET_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 607b77ba244SStefano Zampini 6089566063dSJacob Faibussowitsch PetscCall(MatTranspose(mdata->B, reuse, &mdata->Bt)); 6099566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(mdata->Bt, shell->left_work, NULL)); 6109566063dSJacob Faibussowitsch PetscCall(MatAXPY(D, 1.0, mdata->Bt, str)); 611b77ba244SStefano Zampini } else { 612b77ba244SStefano Zampini MatStructure str = newB ? DIFFERENT_NONZERO_PATTERN : SUBSET_NONZERO_PATTERN; 613b77ba244SStefano Zampini 6149566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(mdata->B, shell->left_work, NULL)); 6159566063dSJacob Faibussowitsch PetscCall(MatAXPY(D, 1.0, mdata->B, str)); 616b77ba244SStefano Zampini } 617b77ba244SStefano Zampini } 618b77ba244SStefano Zampini } 619b77ba244SStefano Zampini break; 620b77ba244SStefano Zampini case MATPRODUCT_AtB: /* s R A^t L B + v R L B + R D L B */ 621b77ba244SStefano Zampini if (shell->right) { 6229566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(D, shell->right, NULL)); 623b77ba244SStefano Zampini if (shell->dshift || shell->vshift != zero) { 624b77ba244SStefano Zampini MatStructure str = newB ? DIFFERENT_NONZERO_PATTERN : SUBSET_NONZERO_PATTERN; 625b77ba244SStefano Zampini 6269566063dSJacob Faibussowitsch if (!shell->right_work) PetscCall(MatCreateVecs(A, &shell->right_work, NULL)); 627b77ba244SStefano Zampini if (shell->dshift) { 6289566063dSJacob Faibussowitsch PetscCall(VecCopy(shell->dshift, shell->right_work)); 6299566063dSJacob Faibussowitsch PetscCall(VecShift(shell->right_work, shell->vshift)); 6309566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->right_work, shell->right_work, shell->right)); 631b77ba244SStefano Zampini } else { 6329566063dSJacob Faibussowitsch PetscCall(VecSet(shell->right_work, shell->vshift)); 633b77ba244SStefano Zampini } 6349566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(mdata->B, shell->right_work, NULL)); 6359566063dSJacob Faibussowitsch PetscCall(MatAXPY(D, 1.0, mdata->B, str)); 636b77ba244SStefano Zampini } 637b77ba244SStefano Zampini } 638b77ba244SStefano Zampini break; 639b77ba244SStefano Zampini case MATPRODUCT_PtAP: /* s B^t L A R B + v B^t L R B + B^t L D R B */ 640b77ba244SStefano Zampini case MATPRODUCT_RARt: /* s B L A R B^t + v B L R B^t + B L D R B^t */ 6419371c9d4SSatish Balay PetscCheck(!shell->dshift && shell->vshift == zero, PetscObjectComm((PetscObject)D), PETSC_ERR_SUP, "MatProductSymbolic type %s not supported for %s and %s matrices with diagonal shift", MatProductTypes[product->type], ((PetscObject)A)->type_name, 6429371c9d4SSatish Balay ((PetscObject)B)->type_name); 643b77ba244SStefano Zampini break; 64498921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)D), PETSC_ERR_SUP, "MatProductSymbolic type %s not supported for %s and %s matrices", MatProductTypes[product->type], ((PetscObject)A)->type_name, ((PetscObject)B)->type_name); 645b77ba244SStefano Zampini } 646b77ba244SStefano Zampini if (shell->axpy && shell->axpy_vscale != zero) { 647b77ba244SStefano Zampini Mat X; 648b77ba244SStefano Zampini PetscObjectState axpy_state; 649b77ba244SStefano Zampini MatStructure str = DIFFERENT_NONZERO_PATTERN; /* not sure it is safe to ever use SUBSET_NONZERO_PATTERN */ 650b77ba244SStefano Zampini 6519566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(shell->axpy, &X)); 6529566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X, &axpy_state)); 65308401ef6SPierre Jolivet PetscCheck(shell->axpy_state == axpy_state, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)"); 654b77ba244SStefano Zampini if (!mdata->axpy) { 655b77ba244SStefano Zampini str = DIFFERENT_NONZERO_PATTERN; 6569566063dSJacob Faibussowitsch PetscCall(MatProductCreate(shell->axpy, B, NULL, &mdata->axpy)); 6579566063dSJacob Faibussowitsch PetscCall(MatProductSetType(mdata->axpy, product->type)); 6589566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(mdata->axpy)); 6599566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(mdata->axpy)); 660b77ba244SStefano Zampini } else { /* May be that shell->axpy has changed */ 661b77ba244SStefano Zampini PetscBool flg; 662b77ba244SStefano Zampini 6639566063dSJacob Faibussowitsch PetscCall(MatProductReplaceMats(shell->axpy, B, NULL, mdata->axpy)); 6649566063dSJacob Faibussowitsch PetscCall(MatHasOperation(mdata->axpy, MATOP_PRODUCTSYMBOLIC, &flg)); 665b77ba244SStefano Zampini if (!flg) { 666b77ba244SStefano Zampini str = DIFFERENT_NONZERO_PATTERN; 6679566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(mdata->axpy)); 6689566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(mdata->axpy)); 669b77ba244SStefano Zampini } 670b77ba244SStefano Zampini } 6719566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(mdata->axpy)); 6729566063dSJacob Faibussowitsch PetscCall(MatAXPY(D, shell->axpy_vscale, mdata->axpy, str)); 673b77ba244SStefano Zampini } 674b77ba244SStefano Zampini } 675b77ba244SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "Missing numeric operation"); 676b77ba244SStefano Zampini PetscFunctionReturn(0); 677b77ba244SStefano Zampini } 678b77ba244SStefano Zampini 6799371c9d4SSatish Balay static PetscErrorCode MatProductSymbolic_Shell_X(Mat D) { 680b77ba244SStefano Zampini Mat_Product *product; 681b77ba244SStefano Zampini Mat A, B; 682b77ba244SStefano Zampini MatShellMatFunctionList matmat; 683b77ba244SStefano Zampini Mat_Shell *shell; 684b77ba244SStefano Zampini PetscBool flg; 685b77ba244SStefano Zampini char composedname[256]; 686b77ba244SStefano Zampini MatMatDataShell *mdata; 687b77ba244SStefano Zampini 688b77ba244SStefano Zampini PetscFunctionBegin; 689b77ba244SStefano Zampini MatCheckProduct(D, 1); 690b77ba244SStefano Zampini product = D->product; 69128b400f6SJacob Faibussowitsch PetscCheck(!product->data, PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "Product data not empty"); 692b77ba244SStefano Zampini A = product->A; 693b77ba244SStefano Zampini B = product->B; 694b77ba244SStefano Zampini shell = (Mat_Shell *)A->data; 695b77ba244SStefano Zampini matmat = shell->matmat; 6969566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(composedname, sizeof(composedname), "MatProductSetFromOptions_%s_%s_C", ((PetscObject)A)->type_name, ((PetscObject)B)->type_name)); 697b77ba244SStefano Zampini while (matmat) { 6989566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(composedname, matmat->composedname, &flg)); 699b77ba244SStefano Zampini flg = (PetscBool)(flg && (matmat->ptype == product->type)); 700b77ba244SStefano Zampini if (flg) break; 701b77ba244SStefano Zampini matmat = matmat->next; 702b77ba244SStefano Zampini } 70328b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Composedname \"%s\" for product type %s not found", composedname, MatProductTypes[product->type]); 704b77ba244SStefano Zampini switch (product->type) { 7059371c9d4SSatish Balay case MATPRODUCT_AB: PetscCall(MatSetSizes(D, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N)); break; 7069371c9d4SSatish Balay case MATPRODUCT_AtB: PetscCall(MatSetSizes(D, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N)); break; 7079371c9d4SSatish Balay case MATPRODUCT_ABt: PetscCall(MatSetSizes(D, A->rmap->n, B->rmap->n, A->rmap->N, B->rmap->N)); break; 7089371c9d4SSatish Balay case MATPRODUCT_RARt: PetscCall(MatSetSizes(D, B->rmap->n, B->rmap->n, B->rmap->N, B->rmap->N)); break; 7099371c9d4SSatish Balay case MATPRODUCT_PtAP: PetscCall(MatSetSizes(D, B->cmap->n, B->cmap->n, B->cmap->N, B->cmap->N)); break; 71098921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)D), PETSC_ERR_SUP, "MatProductSymbolic type %s not supported for %s and %s matrices", MatProductTypes[product->type], ((PetscObject)A)->type_name, ((PetscObject)B)->type_name); 711b77ba244SStefano Zampini } 712b77ba244SStefano Zampini /* respect users who passed in a matrix for which resultname is the base type */ 713b77ba244SStefano Zampini if (matmat->resultname) { 7149566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)D, matmat->resultname, &flg)); 71548a46eb9SPierre Jolivet if (!flg) PetscCall(MatSetType(D, matmat->resultname)); 716b77ba244SStefano Zampini } 717b77ba244SStefano Zampini /* If matrix type was not set or different, we need to reset this pointers */ 718b77ba244SStefano Zampini D->ops->productsymbolic = MatProductSymbolic_Shell_X; 719b77ba244SStefano Zampini D->ops->productnumeric = MatProductNumeric_Shell_X; 720b77ba244SStefano Zampini /* attach product data */ 7219566063dSJacob Faibussowitsch PetscCall(PetscNew(&mdata)); 722b77ba244SStefano Zampini mdata->numeric = matmat->numeric; 723b77ba244SStefano Zampini mdata->destroy = matmat->destroy; 724b77ba244SStefano Zampini if (matmat->symbolic) { 7259566063dSJacob Faibussowitsch PetscCall((*matmat->symbolic)(A, B, D, &mdata->userdata)); 726b77ba244SStefano Zampini } else { /* call general setup if symbolic operation not provided */ 7279566063dSJacob Faibussowitsch PetscCall(MatSetUp(D)); 728b77ba244SStefano Zampini } 72928b400f6SJacob Faibussowitsch PetscCheck(D->product, PetscObjectComm((PetscObject)D), PETSC_ERR_COR, "Product disappeared after user symbolic phase"); 73028b400f6SJacob Faibussowitsch PetscCheck(!D->product->data, PetscObjectComm((PetscObject)D), PETSC_ERR_COR, "Product data not empty after user symbolic phase"); 731b77ba244SStefano Zampini D->product->data = mdata; 732b77ba244SStefano Zampini D->product->destroy = DestroyMatMatDataShell; 733b77ba244SStefano Zampini /* Be sure to reset these pointers if the user did something unexpected */ 734b77ba244SStefano Zampini D->ops->productsymbolic = MatProductSymbolic_Shell_X; 735b77ba244SStefano Zampini D->ops->productnumeric = MatProductNumeric_Shell_X; 736b77ba244SStefano Zampini PetscFunctionReturn(0); 737b77ba244SStefano Zampini } 738b77ba244SStefano Zampini 7399371c9d4SSatish Balay static PetscErrorCode MatProductSetFromOptions_Shell_X(Mat D) { 740b77ba244SStefano Zampini Mat_Product *product; 741b77ba244SStefano Zampini Mat A, B; 742b77ba244SStefano Zampini MatShellMatFunctionList matmat; 743b77ba244SStefano Zampini Mat_Shell *shell; 744b77ba244SStefano Zampini PetscBool flg; 745b77ba244SStefano Zampini char composedname[256]; 746b77ba244SStefano Zampini 747b77ba244SStefano Zampini PetscFunctionBegin; 748b77ba244SStefano Zampini MatCheckProduct(D, 1); 749b77ba244SStefano Zampini product = D->product; 750b77ba244SStefano Zampini A = product->A; 751b77ba244SStefano Zampini B = product->B; 7529566063dSJacob Faibussowitsch PetscCall(MatIsShell(A, &flg)); 753b77ba244SStefano Zampini if (!flg) PetscFunctionReturn(0); 754b77ba244SStefano Zampini shell = (Mat_Shell *)A->data; 755b77ba244SStefano Zampini matmat = shell->matmat; 7569566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(composedname, sizeof(composedname), "MatProductSetFromOptions_%s_%s_C", ((PetscObject)A)->type_name, ((PetscObject)B)->type_name)); 757b77ba244SStefano Zampini while (matmat) { 7589566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(composedname, matmat->composedname, &flg)); 759b77ba244SStefano Zampini flg = (PetscBool)(flg && (matmat->ptype == product->type)); 760b77ba244SStefano Zampini if (flg) break; 761b77ba244SStefano Zampini matmat = matmat->next; 762b77ba244SStefano Zampini } 7639371c9d4SSatish Balay if (flg) { 7649371c9d4SSatish Balay D->ops->productsymbolic = MatProductSymbolic_Shell_X; 7659371c9d4SSatish Balay } else PetscCall(PetscInfo(D, " symbolic product %s not registered for product type %s\n", composedname, MatProductTypes[product->type])); 766b77ba244SStefano Zampini PetscFunctionReturn(0); 767b77ba244SStefano Zampini } 768b77ba244SStefano Zampini 7699371c9d4SSatish Balay static PetscErrorCode MatShellSetMatProductOperation_Private(Mat A, MatProductType ptype, PetscErrorCode (*symbolic)(Mat, Mat, Mat, void **), PetscErrorCode (*numeric)(Mat, Mat, Mat, void *), PetscErrorCode (*destroy)(void *), char *composedname, const char *resultname) { 770b77ba244SStefano Zampini PetscBool flg; 771b77ba244SStefano Zampini Mat_Shell *shell; 772b77ba244SStefano Zampini MatShellMatFunctionList matmat; 773b77ba244SStefano Zampini 774b77ba244SStefano Zampini PetscFunctionBegin; 77528b400f6SJacob Faibussowitsch PetscCheck(numeric, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing numeric routine"); 77628b400f6SJacob Faibussowitsch PetscCheck(composedname, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing composed name"); 777b77ba244SStefano Zampini 778b77ba244SStefano Zampini /* add product callback */ 779b77ba244SStefano Zampini shell = (Mat_Shell *)A->data; 780b77ba244SStefano Zampini matmat = shell->matmat; 781b77ba244SStefano Zampini if (!matmat) { 7829566063dSJacob Faibussowitsch PetscCall(PetscNew(&shell->matmat)); 783b77ba244SStefano Zampini matmat = shell->matmat; 784b77ba244SStefano Zampini } else { 785b77ba244SStefano Zampini MatShellMatFunctionList entry = matmat; 786b77ba244SStefano Zampini while (entry) { 7879566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(composedname, entry->composedname, &flg)); 788b77ba244SStefano Zampini flg = (PetscBool)(flg && (entry->ptype == ptype)); 789b77ba244SStefano Zampini matmat = entry; 7902e956fe4SStefano Zampini if (flg) goto set; 791b77ba244SStefano Zampini entry = entry->next; 792b77ba244SStefano Zampini } 7939566063dSJacob Faibussowitsch PetscCall(PetscNew(&matmat->next)); 794b77ba244SStefano Zampini matmat = matmat->next; 795b77ba244SStefano Zampini } 796b77ba244SStefano Zampini 797843d480fSStefano Zampini set: 798b77ba244SStefano Zampini matmat->symbolic = symbolic; 799b77ba244SStefano Zampini matmat->numeric = numeric; 800b77ba244SStefano Zampini matmat->destroy = destroy; 801b77ba244SStefano Zampini matmat->ptype = ptype; 8029566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat->composedname)); 8039566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat->resultname)); 8049566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(composedname, &matmat->composedname)); 8059566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(resultname, &matmat->resultname)); 8069566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Composing %s for product type %s with result %s\n", matmat->composedname, MatProductTypes[matmat->ptype], matmat->resultname ? matmat->resultname : "not specified")); 8079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, matmat->composedname, MatProductSetFromOptions_Shell_X)); 808b77ba244SStefano Zampini PetscFunctionReturn(0); 809b77ba244SStefano Zampini } 810b77ba244SStefano Zampini 811b77ba244SStefano Zampini /*@C 81211a5261eSBarry Smith MatShellSetMatProductOperation - Allows user to set a matrix matrix operation for a `MATSHELL` shell matrix. 813b77ba244SStefano Zampini 81411a5261eSBarry Smith Logically Collective on A 815b77ba244SStefano Zampini 816b77ba244SStefano Zampini Input Parameters: 81711a5261eSBarry Smith + A - the `MATSHELL` shell matrix 818b77ba244SStefano Zampini . ptype - the product type 819b77ba244SStefano Zampini . symbolic - the function for the symbolic phase (can be NULL) 820b77ba244SStefano Zampini . numeric - the function for the numerical phase 821b77ba244SStefano Zampini . destroy - the function for the destruction of the needed data generated during the symbolic phase (can be NULL) 822b77ba244SStefano Zampini . Btype - the matrix type for the matrix to be multiplied against 823b77ba244SStefano Zampini - Ctype - the matrix type for the result (can be NULL) 824b77ba244SStefano Zampini 825b77ba244SStefano Zampini Level: advanced 826b77ba244SStefano Zampini 827b77ba244SStefano Zampini Usage: 828b77ba244SStefano Zampini $ extern PetscErrorCode usersymbolic(Mat,Mat,Mat,void**); 829b77ba244SStefano Zampini $ extern PetscErrorCode usernumeric(Mat,Mat,Mat,void*); 830b77ba244SStefano Zampini $ extern PetscErrorCode userdestroy(void*); 831b77ba244SStefano Zampini $ MatCreateShell(comm,m,n,M,N,ctx,&A); 832b77ba244SStefano Zampini $ MatShellSetMatProductOperation(A,MATPRODUCT_AB,usersymbolic,usernumeric,userdestroy,MATSEQAIJ,MATDENSE); 833b77ba244SStefano Zampini $ [ create B of type SEQAIJ etc..] 834b77ba244SStefano Zampini $ MatProductCreate(A,B,NULL,&C); 835b77ba244SStefano Zampini $ MatProductSetType(C,MATPRODUCT_AB); 836b77ba244SStefano Zampini $ MatProductSetFromOptions(C); 837b77ba244SStefano Zampini $ MatProductSymbolic(C); -> actually runs the user defined symbolic operation 838b77ba244SStefano Zampini $ MatProductNumeric(C); -> actually runs the user defined numeric operation 839b77ba244SStefano Zampini $ [ use C = A*B ] 840b77ba244SStefano Zampini 841b77ba244SStefano Zampini Notes: 84211a5261eSBarry Smith `MATPRODUCT_ABC` is not supported yet. Not supported in Fortran. 84311a5261eSBarry Smith 84411a5261eSBarry Smith If the symbolic phase is not specified, `MatSetUp()` is called on the result matrix that must have its type set if Ctype is NULL. 84511a5261eSBarry Smith 846b77ba244SStefano Zampini Any additional data needed by the matrix product needs to be returned during the symbolic phase and destroyed with the destroy callback. 847b77ba244SStefano Zampini PETSc will take care of calling the user-defined callbacks. 848b77ba244SStefano Zampini It is allowed to specify the same callbacks for different Btype matrix types. 849b77ba244SStefano Zampini The couple (Btype,ptype) uniquely identifies the operation: the last specified callbacks takes precedence. 850b77ba244SStefano Zampini 85111a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellSetContext()`, `MatSetOperation()`, `MatProductType`, `MatType`, `MatSetUp()` 852b77ba244SStefano Zampini @*/ 8539371c9d4SSatish Balay PetscErrorCode MatShellSetMatProductOperation(Mat A, MatProductType ptype, PetscErrorCode (*symbolic)(Mat, Mat, Mat, void **), PetscErrorCode (*numeric)(Mat, Mat, Mat, void *), PetscErrorCode (*destroy)(void *), MatType Btype, MatType Ctype) { 854b77ba244SStefano Zampini PetscFunctionBegin; 855b77ba244SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 856b77ba244SStefano Zampini PetscValidLogicalCollectiveEnum(A, ptype, 2); 85708401ef6SPierre Jolivet PetscCheck(ptype != MATPRODUCT_ABC, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for product type %s", MatProductTypes[ptype]); 85828b400f6SJacob Faibussowitsch PetscCheck(numeric, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Missing numeric routine, argument 4"); 859dadcf809SJacob Faibussowitsch PetscValidCharPointer(Btype, 6); 860dadcf809SJacob Faibussowitsch if (Ctype) PetscValidCharPointer(Ctype, 7); 861cac4c232SBarry Smith PetscTryMethod(A, "MatShellSetMatProductOperation_C", (Mat, MatProductType, PetscErrorCode(*)(Mat, Mat, Mat, void **), PetscErrorCode(*)(Mat, Mat, Mat, void *), PetscErrorCode(*)(void *), MatType, MatType), (A, ptype, symbolic, numeric, destroy, Btype, Ctype)); 862b77ba244SStefano Zampini PetscFunctionReturn(0); 863b77ba244SStefano Zampini } 864b77ba244SStefano Zampini 8659371c9d4SSatish Balay static PetscErrorCode MatShellSetMatProductOperation_Shell(Mat A, MatProductType ptype, PetscErrorCode (*symbolic)(Mat, Mat, Mat, void **), PetscErrorCode (*numeric)(Mat, Mat, Mat, void *), PetscErrorCode (*destroy)(void *), MatType Btype, MatType Ctype) { 866b77ba244SStefano Zampini PetscBool flg; 867b77ba244SStefano Zampini char composedname[256]; 868b77ba244SStefano Zampini MatRootName Bnames = MatRootNameList, Cnames = MatRootNameList; 869b77ba244SStefano Zampini PetscMPIInt size; 870b77ba244SStefano Zampini 871b77ba244SStefano Zampini PetscFunctionBegin; 872b77ba244SStefano Zampini PetscValidType(A, 1); 873b77ba244SStefano Zampini while (Bnames) { /* user passed in the root name */ 8749566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(Btype, Bnames->rname, &flg)); 875b77ba244SStefano Zampini if (flg) break; 876b77ba244SStefano Zampini Bnames = Bnames->next; 877b77ba244SStefano Zampini } 878b77ba244SStefano Zampini while (Cnames) { /* user passed in the root name */ 8799566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(Ctype, Cnames->rname, &flg)); 880b77ba244SStefano Zampini if (flg) break; 881b77ba244SStefano Zampini Cnames = Cnames->next; 882b77ba244SStefano Zampini } 8839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 884b77ba244SStefano Zampini Btype = Bnames ? (size > 1 ? Bnames->mname : Bnames->sname) : Btype; 885b77ba244SStefano Zampini Ctype = Cnames ? (size > 1 ? Cnames->mname : Cnames->sname) : Ctype; 8869566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(composedname, sizeof(composedname), "MatProductSetFromOptions_%s_%s_C", ((PetscObject)A)->type_name, Btype)); 8879566063dSJacob Faibussowitsch PetscCall(MatShellSetMatProductOperation_Private(A, ptype, symbolic, numeric, destroy, composedname, Ctype)); 8883a40ed3dSBarry Smith PetscFunctionReturn(0); 889e51e0e81SBarry Smith } 890e51e0e81SBarry Smith 8919371c9d4SSatish Balay static PetscErrorCode MatCopy_Shell(Mat A, Mat B, MatStructure str) { 89228f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell *)A->data, *shellB = (Mat_Shell *)B->data; 893cb8c8a77SBarry Smith PetscBool matflg; 894b77ba244SStefano Zampini MatShellMatFunctionList matmatA; 8957fabda3fSAlex Fikl 8967fabda3fSAlex Fikl PetscFunctionBegin; 8979566063dSJacob Faibussowitsch PetscCall(MatIsShell(B, &matflg)); 89828b400f6SJacob Faibussowitsch PetscCheck(matflg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix %s not derived from MATSHELL", ((PetscObject)B)->type_name); 8997fabda3fSAlex Fikl 9009566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(B->ops, A->ops, sizeof(struct _MatOps))); 9019566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(shellB->ops, shellA->ops, sizeof(struct _MatShellOps))); 9027fabda3fSAlex Fikl 9031baa6e33SBarry Smith if (shellA->ops->copy) PetscCall((*shellA->ops->copy)(A, B, str)); 9047fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 9057fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 9060c0fd78eSBarry Smith if (shellA->dshift) { 90748a46eb9SPierre Jolivet if (!shellB->dshift) PetscCall(VecDuplicate(shellA->dshift, &shellB->dshift)); 9089566063dSJacob Faibussowitsch PetscCall(VecCopy(shellA->dshift, shellB->dshift)); 9097fabda3fSAlex Fikl } else { 9109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shellB->dshift)); 9117fabda3fSAlex Fikl } 9120c0fd78eSBarry Smith if (shellA->left) { 91348a46eb9SPierre Jolivet if (!shellB->left) PetscCall(VecDuplicate(shellA->left, &shellB->left)); 9149566063dSJacob Faibussowitsch PetscCall(VecCopy(shellA->left, shellB->left)); 9157fabda3fSAlex Fikl } else { 9169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shellB->left)); 9177fabda3fSAlex Fikl } 9180c0fd78eSBarry Smith if (shellA->right) { 91948a46eb9SPierre Jolivet if (!shellB->right) PetscCall(VecDuplicate(shellA->right, &shellB->right)); 9209566063dSJacob Faibussowitsch PetscCall(VecCopy(shellA->right, shellB->right)); 9217fabda3fSAlex Fikl } else { 9229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shellB->right)); 9237fabda3fSAlex Fikl } 9249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&shellB->axpy)); 925ef5c7bd2SStefano Zampini shellB->axpy_vscale = 0.0; 926ef5c7bd2SStefano Zampini shellB->axpy_state = 0; 92740e381c3SBarry Smith if (shellA->axpy) { 9289566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)shellA->axpy)); 92940e381c3SBarry Smith shellB->axpy = shellA->axpy; 93040e381c3SBarry Smith shellB->axpy_vscale = shellA->axpy_vscale; 931ef5c7bd2SStefano Zampini shellB->axpy_state = shellA->axpy_state; 93240e381c3SBarry Smith } 93345960306SStefano Zampini if (shellA->zrows) { 9349566063dSJacob Faibussowitsch PetscCall(ISDuplicate(shellA->zrows, &shellB->zrows)); 93548a46eb9SPierre Jolivet if (shellA->zcols) PetscCall(ISDuplicate(shellA->zcols, &shellB->zcols)); 9369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shellA->zvals, &shellB->zvals)); 9379566063dSJacob Faibussowitsch PetscCall(VecCopy(shellA->zvals, shellB->zvals)); 9389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shellA->zvals_w, &shellB->zvals_w)); 9399566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)shellA->zvals_sct_r)); 9409566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)shellA->zvals_sct_c)); 94145960306SStefano Zampini shellB->zvals_sct_r = shellA->zvals_sct_r; 94245960306SStefano Zampini shellB->zvals_sct_c = shellA->zvals_sct_c; 94345960306SStefano Zampini } 944b77ba244SStefano Zampini 945b77ba244SStefano Zampini matmatA = shellA->matmat; 946b77ba244SStefano Zampini if (matmatA) { 947b77ba244SStefano Zampini while (matmatA->next) { 9489566063dSJacob Faibussowitsch PetscCall(MatShellSetMatProductOperation_Private(B, matmatA->ptype, matmatA->symbolic, matmatA->numeric, matmatA->destroy, matmatA->composedname, matmatA->resultname)); 949b77ba244SStefano Zampini matmatA = matmatA->next; 950b77ba244SStefano Zampini } 951b77ba244SStefano Zampini } 9527fabda3fSAlex Fikl PetscFunctionReturn(0); 9537fabda3fSAlex Fikl } 9547fabda3fSAlex Fikl 9559371c9d4SSatish Balay static PetscErrorCode MatDuplicate_Shell(Mat mat, MatDuplicateOption op, Mat *M) { 956cb8c8a77SBarry Smith PetscFunctionBegin; 957800f99ffSJeremy L Thompson PetscCall(MatCreateShell(PetscObjectComm((PetscObject)mat), mat->rmap->n, mat->cmap->n, mat->rmap->N, mat->cmap->N, NULL, M)); 958800f99ffSJeremy L Thompson ((Mat_Shell *)(*M)->data)->ctxcontainer = ((Mat_Shell *)mat->data)->ctxcontainer; 959800f99ffSJeremy L Thompson PetscCall(PetscObjectCompose((PetscObject)(*M), "MatShell ctx", (PetscObject)((Mat_Shell *)(*M)->data)->ctxcontainer)); 9609566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)(*M), ((PetscObject)mat)->type_name)); 96148a46eb9SPierre Jolivet if (op != MAT_DO_NOT_COPY_VALUES) PetscCall(MatCopy(mat, *M, SAME_NONZERO_PATTERN)); 962cb8c8a77SBarry Smith PetscFunctionReturn(0); 963cb8c8a77SBarry Smith } 964cb8c8a77SBarry Smith 9659371c9d4SSatish Balay PetscErrorCode MatMult_Shell(Mat A, Vec x, Vec y) { 966ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell *)A->data; 96725578ef6SJed Brown Vec xx; 968e3f487b0SBarry Smith PetscObjectState instate, outstate; 969ef66eb69SBarry Smith 970ef66eb69SBarry Smith PetscFunctionBegin; 9719566063dSJacob Faibussowitsch PetscCall(MatShellPreZeroRight(A, x, &xx)); 9729566063dSJacob Faibussowitsch PetscCall(MatShellPreScaleRight(A, xx, &xx)); 9739566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 9749566063dSJacob Faibussowitsch PetscCall((*shell->ops->mult)(A, xx, y)); 9759566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 976e3f487b0SBarry Smith if (instate == outstate) { 977e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 9789566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)y)); 979e3f487b0SBarry Smith } 9809566063dSJacob Faibussowitsch PetscCall(MatShellShiftAndScale(A, xx, y)); 9819566063dSJacob Faibussowitsch PetscCall(MatShellPostScaleLeft(A, y)); 9829566063dSJacob Faibussowitsch PetscCall(MatShellPostZeroLeft(A, y)); 9839f137db4SBarry Smith 9849f137db4SBarry Smith if (shell->axpy) { 985ef5c7bd2SStefano Zampini Mat X; 986ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 987ef5c7bd2SStefano Zampini 9889566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(shell->axpy, &X)); 9899566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X, &axpy_state)); 99008401ef6SPierre Jolivet PetscCheck(shell->axpy_state == axpy_state, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)"); 991b77ba244SStefano Zampini 9929566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(shell->axpy, shell->axpy_right ? NULL : &shell->axpy_right, shell->axpy_left ? NULL : &shell->axpy_left)); 9939566063dSJacob Faibussowitsch PetscCall(VecCopy(x, shell->axpy_right)); 9949566063dSJacob Faibussowitsch PetscCall(MatMult(shell->axpy, shell->axpy_right, shell->axpy_left)); 9959566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, shell->axpy_vscale, shell->axpy_left)); 9969f137db4SBarry Smith } 997ef66eb69SBarry Smith PetscFunctionReturn(0); 998ef66eb69SBarry Smith } 999ef66eb69SBarry Smith 10009371c9d4SSatish Balay static PetscErrorCode MatMultAdd_Shell(Mat A, Vec x, Vec y, Vec z) { 10015edf6516SJed Brown Mat_Shell *shell = (Mat_Shell *)A->data; 10025edf6516SJed Brown 10035edf6516SJed Brown PetscFunctionBegin; 10045edf6516SJed Brown if (y == z) { 10059566063dSJacob Faibussowitsch if (!shell->right_add_work) PetscCall(VecDuplicate(z, &shell->right_add_work)); 10069566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, shell->right_add_work)); 10079566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, 1.0, shell->right_add_work)); 10085edf6516SJed Brown } else { 10099566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, z)); 10109566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, 1.0, y)); 10115edf6516SJed Brown } 10125edf6516SJed Brown PetscFunctionReturn(0); 10135edf6516SJed Brown } 10145edf6516SJed Brown 10159371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_Shell(Mat A, Vec x, Vec y) { 101618be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell *)A->data; 101725578ef6SJed Brown Vec xx; 1018e3f487b0SBarry Smith PetscObjectState instate, outstate; 101918be62a5SSatish Balay 102018be62a5SSatish Balay PetscFunctionBegin; 10219566063dSJacob Faibussowitsch PetscCall(MatShellPreZeroLeft(A, x, &xx)); 10229566063dSJacob Faibussowitsch PetscCall(MatShellPreScaleLeft(A, xx, &xx)); 10239566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 10249566063dSJacob Faibussowitsch PetscCall((*shell->ops->multtranspose)(A, xx, y)); 10259566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 1026e3f487b0SBarry Smith if (instate == outstate) { 1027e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 10289566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)y)); 1029e3f487b0SBarry Smith } 10309566063dSJacob Faibussowitsch PetscCall(MatShellShiftAndScale(A, xx, y)); 10319566063dSJacob Faibussowitsch PetscCall(MatShellPostScaleRight(A, y)); 10329566063dSJacob Faibussowitsch PetscCall(MatShellPostZeroRight(A, y)); 1033350ec83bSStefano Zampini 1034350ec83bSStefano Zampini if (shell->axpy) { 1035ef5c7bd2SStefano Zampini Mat X; 1036ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 1037ef5c7bd2SStefano Zampini 10389566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(shell->axpy, &X)); 10399566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X, &axpy_state)); 104008401ef6SPierre Jolivet PetscCheck(shell->axpy_state == axpy_state, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)"); 10419566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(shell->axpy, shell->axpy_right ? NULL : &shell->axpy_right, shell->axpy_left ? NULL : &shell->axpy_left)); 10429566063dSJacob Faibussowitsch PetscCall(VecCopy(x, shell->axpy_left)); 10439566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(shell->axpy, shell->axpy_left, shell->axpy_right)); 10449566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, shell->axpy_vscale, shell->axpy_right)); 1045350ec83bSStefano Zampini } 104618be62a5SSatish Balay PetscFunctionReturn(0); 104718be62a5SSatish Balay } 104818be62a5SSatish Balay 10499371c9d4SSatish Balay static PetscErrorCode MatMultTransposeAdd_Shell(Mat A, Vec x, Vec y, Vec z) { 10505edf6516SJed Brown Mat_Shell *shell = (Mat_Shell *)A->data; 10515edf6516SJed Brown 10525edf6516SJed Brown PetscFunctionBegin; 10535edf6516SJed Brown if (y == z) { 10549566063dSJacob Faibussowitsch if (!shell->left_add_work) PetscCall(VecDuplicate(z, &shell->left_add_work)); 10559566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, shell->left_add_work)); 10569566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, 1.0, shell->left_add_work)); 10575edf6516SJed Brown } else { 10589566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, z)); 10599566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, 1.0, y)); 10605edf6516SJed Brown } 10615edf6516SJed Brown PetscFunctionReturn(0); 10625edf6516SJed Brown } 10635edf6516SJed Brown 10640c0fd78eSBarry Smith /* 10650c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 10660c0fd78eSBarry Smith */ 10679371c9d4SSatish Balay static PetscErrorCode MatGetDiagonal_Shell(Mat A, Vec v) { 106818be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell *)A->data; 106918be62a5SSatish Balay 107018be62a5SSatish Balay PetscFunctionBegin; 10710c0fd78eSBarry Smith if (shell->ops->getdiagonal) { 10729566063dSJacob Faibussowitsch PetscCall((*shell->ops->getdiagonal)(A, v)); 1073305211d5SBarry 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,...)"); 10749566063dSJacob Faibussowitsch PetscCall(VecScale(v, shell->vscale)); 10751baa6e33SBarry Smith if (shell->dshift) PetscCall(VecAXPY(v, 1.0, shell->dshift)); 10769566063dSJacob Faibussowitsch PetscCall(VecShift(v, shell->vshift)); 10779566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseMult(v, v, shell->left)); 10789566063dSJacob Faibussowitsch if (shell->right) PetscCall(VecPointwiseMult(v, v, shell->right)); 107945960306SStefano Zampini if (shell->zrows) { 10809566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_r, shell->zvals, v, INSERT_VALUES, SCATTER_REVERSE)); 10819566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_r, shell->zvals, v, INSERT_VALUES, SCATTER_REVERSE)); 108245960306SStefano Zampini } 108381c02519SBarry Smith if (shell->axpy) { 1084ef5c7bd2SStefano Zampini Mat X; 1085ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 1086ef5c7bd2SStefano Zampini 10879566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(shell->axpy, &X)); 10889566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X, &axpy_state)); 108908401ef6SPierre Jolivet PetscCheck(shell->axpy_state == axpy_state, PetscObjectComm((PetscObject)A), PETSC_ERR_ORDER, "Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)"); 10909566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(shell->axpy, NULL, shell->axpy_left ? NULL : &shell->axpy_left)); 10919566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(shell->axpy, shell->axpy_left)); 10929566063dSJacob Faibussowitsch PetscCall(VecAXPY(v, shell->axpy_vscale, shell->axpy_left)); 109381c02519SBarry Smith } 109418be62a5SSatish Balay PetscFunctionReturn(0); 109518be62a5SSatish Balay } 109618be62a5SSatish Balay 10979371c9d4SSatish Balay static PetscErrorCode MatShift_Shell(Mat Y, PetscScalar a) { 1098ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell *)Y->data; 1099789d8953SBarry Smith PetscBool flg; 1100b24ad042SBarry Smith 1101ef66eb69SBarry Smith PetscFunctionBegin; 11029566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(Y, &flg)); 110328b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)Y), PETSC_ERR_SUP, "Cannot shift shell matrix if it is not congruent"); 11040c0fd78eSBarry Smith if (shell->left || shell->right) { 11058fe8eb27SJed Brown if (!shell->dshift) { 11069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift)); 11079566063dSJacob Faibussowitsch PetscCall(VecSet(shell->dshift, a)); 11080c0fd78eSBarry Smith } else { 11099566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseMult(shell->dshift, shell->dshift, shell->left)); 11109566063dSJacob Faibussowitsch if (shell->right) PetscCall(VecPointwiseMult(shell->dshift, shell->dshift, shell->right)); 11119566063dSJacob Faibussowitsch PetscCall(VecShift(shell->dshift, a)); 11120c0fd78eSBarry Smith } 11139566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseDivide(shell->dshift, shell->dshift, shell->left)); 11149566063dSJacob Faibussowitsch if (shell->right) PetscCall(VecPointwiseDivide(shell->dshift, shell->dshift, shell->right)); 11158fe8eb27SJed Brown } else shell->vshift += a; 11161baa6e33SBarry Smith if (shell->zrows) PetscCall(VecShift(shell->zvals, a)); 1117ef66eb69SBarry Smith PetscFunctionReturn(0); 1118ef66eb69SBarry Smith } 1119ef66eb69SBarry Smith 11209371c9d4SSatish Balay static PetscErrorCode MatDiagonalSet_Shell_Private(Mat A, Vec D, PetscScalar s) { 11216464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell *)A->data; 11226464bf51SAlex Fikl 11236464bf51SAlex Fikl PetscFunctionBegin; 11249566063dSJacob Faibussowitsch if (!shell->dshift) PetscCall(VecDuplicate(D, &shell->dshift)); 11250c0fd78eSBarry Smith if (shell->left || shell->right) { 11269566063dSJacob Faibussowitsch if (!shell->right_work) PetscCall(VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work)); 11279f137db4SBarry Smith if (shell->left && shell->right) { 11289566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(shell->right_work, D, shell->left)); 11299566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(shell->right_work, shell->right_work, shell->right)); 11309f137db4SBarry Smith } else if (shell->left) { 11319566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(shell->right_work, D, shell->left)); 11329f137db4SBarry Smith } else { 11339566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(shell->right_work, D, shell->right)); 11349f137db4SBarry Smith } 11359566063dSJacob Faibussowitsch PetscCall(VecAXPY(shell->dshift, s, shell->right_work)); 11360c0fd78eSBarry Smith } else { 11379566063dSJacob Faibussowitsch PetscCall(VecAXPY(shell->dshift, s, D)); 1138b253ae0bS“Marek } 1139b253ae0bS“Marek PetscFunctionReturn(0); 1140b253ae0bS“Marek } 1141b253ae0bS“Marek 11429371c9d4SSatish Balay PetscErrorCode MatDiagonalSet_Shell(Mat A, Vec D, InsertMode ins) { 114345960306SStefano Zampini Mat_Shell *shell = (Mat_Shell *)A->data; 1144b253ae0bS“Marek Vec d; 1145789d8953SBarry Smith PetscBool flg; 1146b253ae0bS“Marek 1147b253ae0bS“Marek PetscFunctionBegin; 11489566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A, &flg)); 114928b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot diagonal set or shift shell matrix if it is not congruent"); 1150b253ae0bS“Marek if (ins == INSERT_VALUES) { 11519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(D, &d)); 11529566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, d)); 11539566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet_Shell_Private(A, d, -1.)); 11549566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet_Shell_Private(A, D, 1.)); 11559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d)); 11561baa6e33SBarry Smith if (shell->zrows) PetscCall(VecCopy(D, shell->zvals)); 1157b253ae0bS“Marek } else { 11589566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet_Shell_Private(A, D, 1.)); 11591baa6e33SBarry Smith if (shell->zrows) PetscCall(VecAXPY(shell->zvals, 1.0, D)); 11606464bf51SAlex Fikl } 11616464bf51SAlex Fikl PetscFunctionReturn(0); 11626464bf51SAlex Fikl } 11636464bf51SAlex Fikl 11649371c9d4SSatish Balay static PetscErrorCode MatScale_Shell(Mat Y, PetscScalar a) { 1165ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell *)Y->data; 1166b24ad042SBarry Smith 1167ef66eb69SBarry Smith PetscFunctionBegin; 1168f4df32b1SMatthew Knepley shell->vscale *= a; 11690c0fd78eSBarry Smith shell->vshift *= a; 11701baa6e33SBarry Smith if (shell->dshift) PetscCall(VecScale(shell->dshift, a)); 117181c02519SBarry Smith shell->axpy_vscale *= a; 11721baa6e33SBarry Smith if (shell->zrows) PetscCall(VecScale(shell->zvals, a)); 11738fe8eb27SJed Brown PetscFunctionReturn(0); 117418be62a5SSatish Balay } 11758fe8eb27SJed Brown 11769371c9d4SSatish Balay static PetscErrorCode MatDiagonalScale_Shell(Mat Y, Vec left, Vec right) { 11778fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell *)Y->data; 11788fe8eb27SJed Brown 11798fe8eb27SJed Brown PetscFunctionBegin; 11808fe8eb27SJed Brown if (left) { 11810c0fd78eSBarry Smith if (!shell->left) { 11829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left, &shell->left)); 11839566063dSJacob Faibussowitsch PetscCall(VecCopy(left, shell->left)); 11840c0fd78eSBarry Smith } else { 11859566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->left, shell->left, left)); 118618be62a5SSatish Balay } 11871baa6e33SBarry Smith if (shell->zrows) PetscCall(VecPointwiseMult(shell->zvals, shell->zvals, left)); 1188ef66eb69SBarry Smith } 11898fe8eb27SJed Brown if (right) { 11900c0fd78eSBarry Smith if (!shell->right) { 11919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right, &shell->right)); 11929566063dSJacob Faibussowitsch PetscCall(VecCopy(right, shell->right)); 11930c0fd78eSBarry Smith } else { 11949566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->right, shell->right, right)); 11958fe8eb27SJed Brown } 119645960306SStefano Zampini if (shell->zrows) { 119748a46eb9SPierre Jolivet if (!shell->left_work) PetscCall(MatCreateVecs(Y, NULL, &shell->left_work)); 11989566063dSJacob Faibussowitsch PetscCall(VecSet(shell->zvals_w, 1.0)); 11999566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_c, right, shell->zvals_w, INSERT_VALUES, SCATTER_FORWARD)); 12009566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_c, right, shell->zvals_w, INSERT_VALUES, SCATTER_FORWARD)); 12019566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->zvals, shell->zvals, shell->zvals_w)); 120245960306SStefano Zampini } 12038fe8eb27SJed Brown } 12041baa6e33SBarry Smith if (shell->axpy) PetscCall(MatDiagonalScale(shell->axpy, left, right)); 1205ef66eb69SBarry Smith PetscFunctionReturn(0); 1206ef66eb69SBarry Smith } 1207ef66eb69SBarry Smith 12089371c9d4SSatish Balay PetscErrorCode MatAssemblyEnd_Shell(Mat Y, MatAssemblyType t) { 1209ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell *)Y->data; 1210ef66eb69SBarry Smith 1211ef66eb69SBarry Smith PetscFunctionBegin; 12128fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 1213ef66eb69SBarry Smith shell->vshift = 0.0; 1214ef66eb69SBarry Smith shell->vscale = 1.0; 1215ef5c7bd2SStefano Zampini shell->axpy_vscale = 0.0; 1216ef5c7bd2SStefano Zampini shell->axpy_state = 0; 12179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->dshift)); 12189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->left)); 12199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->right)); 12209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&shell->axpy)); 12219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->axpy_left)); 12229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->axpy_right)); 12239566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_c)); 12249566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_r)); 12259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zrows)); 12269566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zcols)); 1227ef66eb69SBarry Smith } 1228ef66eb69SBarry Smith PetscFunctionReturn(0); 1229ef66eb69SBarry Smith } 1230ef66eb69SBarry Smith 12319371c9d4SSatish Balay static PetscErrorCode MatMissingDiagonal_Shell(Mat A, PetscBool *missing, PetscInt *d) { 12323b49f96aSBarry Smith PetscFunctionBegin; 12333b49f96aSBarry Smith *missing = PETSC_FALSE; 12343b49f96aSBarry Smith PetscFunctionReturn(0); 12353b49f96aSBarry Smith } 12363b49f96aSBarry Smith 12379371c9d4SSatish Balay PetscErrorCode MatAXPY_Shell(Mat Y, PetscScalar a, Mat X, MatStructure str) { 12389f137db4SBarry Smith Mat_Shell *shell = (Mat_Shell *)Y->data; 12399f137db4SBarry Smith 12409f137db4SBarry Smith PetscFunctionBegin; 1241ef5c7bd2SStefano Zampini if (X == Y) { 12429566063dSJacob Faibussowitsch PetscCall(MatScale(Y, 1.0 + a)); 1243ef5c7bd2SStefano Zampini PetscFunctionReturn(0); 1244ef5c7bd2SStefano Zampini } 1245ef5c7bd2SStefano Zampini if (!shell->axpy) { 12469566063dSJacob Faibussowitsch PetscCall(MatConvertFrom_Shell(X, MATSHELL, MAT_INITIAL_MATRIX, &shell->axpy)); 12479f137db4SBarry Smith shell->axpy_vscale = a; 12489566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X, &shell->axpy_state)); 1249ef5c7bd2SStefano Zampini } else { 12509566063dSJacob Faibussowitsch PetscCall(MatAXPY(shell->axpy, a / shell->axpy_vscale, X, str)); 1251ef5c7bd2SStefano Zampini } 12529f137db4SBarry Smith PetscFunctionReturn(0); 12539f137db4SBarry Smith } 12549f137db4SBarry Smith 1255f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL, 1256f4259b30SLisandro Dalcin NULL, 1257f4259b30SLisandro Dalcin NULL, 1258f4259b30SLisandro Dalcin NULL, 12590c0fd78eSBarry Smith /* 4*/ MatMultAdd_Shell, 1260f4259b30SLisandro Dalcin NULL, 12610c0fd78eSBarry Smith MatMultTransposeAdd_Shell, 1262f4259b30SLisandro Dalcin NULL, 1263f4259b30SLisandro Dalcin NULL, 1264f4259b30SLisandro Dalcin NULL, 1265f4259b30SLisandro Dalcin /*10*/ NULL, 1266f4259b30SLisandro Dalcin NULL, 1267f4259b30SLisandro Dalcin NULL, 1268f4259b30SLisandro Dalcin NULL, 1269f4259b30SLisandro Dalcin NULL, 1270f4259b30SLisandro Dalcin /*15*/ NULL, 1271f4259b30SLisandro Dalcin NULL, 1272f4259b30SLisandro Dalcin NULL, 12738fe8eb27SJed Brown MatDiagonalScale_Shell, 1274f4259b30SLisandro Dalcin NULL, 1275f4259b30SLisandro Dalcin /*20*/ NULL, 1276ef66eb69SBarry Smith MatAssemblyEnd_Shell, 1277f4259b30SLisandro Dalcin NULL, 1278f4259b30SLisandro Dalcin NULL, 127945960306SStefano Zampini /*24*/ MatZeroRows_Shell, 1280f4259b30SLisandro Dalcin NULL, 1281f4259b30SLisandro Dalcin NULL, 1282f4259b30SLisandro Dalcin NULL, 1283f4259b30SLisandro Dalcin NULL, 1284f4259b30SLisandro Dalcin /*29*/ NULL, 1285f4259b30SLisandro Dalcin NULL, 1286f4259b30SLisandro Dalcin NULL, 1287f4259b30SLisandro Dalcin NULL, 1288f4259b30SLisandro Dalcin NULL, 1289cb8c8a77SBarry Smith /*34*/ MatDuplicate_Shell, 1290f4259b30SLisandro Dalcin NULL, 1291f4259b30SLisandro Dalcin NULL, 1292f4259b30SLisandro Dalcin NULL, 1293f4259b30SLisandro Dalcin NULL, 12949f137db4SBarry Smith /*39*/ MatAXPY_Shell, 1295f4259b30SLisandro Dalcin NULL, 1296f4259b30SLisandro Dalcin NULL, 1297f4259b30SLisandro Dalcin NULL, 1298cb8c8a77SBarry Smith MatCopy_Shell, 1299f4259b30SLisandro Dalcin /*44*/ NULL, 1300ef66eb69SBarry Smith MatScale_Shell, 1301ef66eb69SBarry Smith MatShift_Shell, 13026464bf51SAlex Fikl MatDiagonalSet_Shell, 130345960306SStefano Zampini MatZeroRowsColumns_Shell, 1304f4259b30SLisandro Dalcin /*49*/ NULL, 1305f4259b30SLisandro Dalcin NULL, 1306f4259b30SLisandro Dalcin NULL, 1307f4259b30SLisandro Dalcin NULL, 1308f4259b30SLisandro Dalcin NULL, 1309f4259b30SLisandro Dalcin /*54*/ NULL, 1310f4259b30SLisandro Dalcin NULL, 1311f4259b30SLisandro Dalcin NULL, 1312f4259b30SLisandro Dalcin NULL, 1313f4259b30SLisandro Dalcin NULL, 1314f4259b30SLisandro Dalcin /*59*/ NULL, 1315b9b97703SBarry Smith MatDestroy_Shell, 1316f4259b30SLisandro Dalcin NULL, 1317251fad3fSStefano Zampini MatConvertFrom_Shell, 1318f4259b30SLisandro Dalcin NULL, 1319f4259b30SLisandro Dalcin /*64*/ NULL, 1320f4259b30SLisandro Dalcin NULL, 1321f4259b30SLisandro Dalcin NULL, 1322f4259b30SLisandro Dalcin NULL, 1323f4259b30SLisandro Dalcin NULL, 1324f4259b30SLisandro Dalcin /*69*/ NULL, 1325f4259b30SLisandro Dalcin NULL, 1326c87e5d42SMatthew Knepley MatConvert_Shell, 1327f4259b30SLisandro Dalcin NULL, 1328f4259b30SLisandro Dalcin NULL, 1329f4259b30SLisandro Dalcin /*74*/ NULL, 1330f4259b30SLisandro Dalcin NULL, 1331f4259b30SLisandro Dalcin NULL, 1332f4259b30SLisandro Dalcin NULL, 1333f4259b30SLisandro Dalcin NULL, 1334f4259b30SLisandro Dalcin /*79*/ NULL, 1335f4259b30SLisandro Dalcin NULL, 1336f4259b30SLisandro Dalcin NULL, 1337f4259b30SLisandro Dalcin NULL, 1338f4259b30SLisandro Dalcin NULL, 1339f4259b30SLisandro Dalcin /*84*/ NULL, 1340f4259b30SLisandro Dalcin NULL, 1341f4259b30SLisandro Dalcin NULL, 1342f4259b30SLisandro Dalcin NULL, 1343f4259b30SLisandro Dalcin NULL, 1344f4259b30SLisandro Dalcin /*89*/ NULL, 1345f4259b30SLisandro Dalcin NULL, 1346f4259b30SLisandro Dalcin NULL, 1347f4259b30SLisandro Dalcin NULL, 1348f4259b30SLisandro Dalcin NULL, 1349f4259b30SLisandro Dalcin /*94*/ NULL, 1350f4259b30SLisandro Dalcin NULL, 1351f4259b30SLisandro Dalcin NULL, 1352f4259b30SLisandro Dalcin NULL, 1353f4259b30SLisandro Dalcin NULL, 1354f4259b30SLisandro Dalcin /*99*/ NULL, 1355f4259b30SLisandro Dalcin NULL, 1356f4259b30SLisandro Dalcin NULL, 1357f4259b30SLisandro Dalcin NULL, 1358f4259b30SLisandro Dalcin NULL, 1359f4259b30SLisandro Dalcin /*104*/ NULL, 1360f4259b30SLisandro Dalcin NULL, 1361f4259b30SLisandro Dalcin NULL, 1362f4259b30SLisandro Dalcin NULL, 1363f4259b30SLisandro Dalcin NULL, 1364f4259b30SLisandro Dalcin /*109*/ NULL, 1365f4259b30SLisandro Dalcin NULL, 1366f4259b30SLisandro Dalcin NULL, 1367f4259b30SLisandro Dalcin NULL, 13683b49f96aSBarry Smith MatMissingDiagonal_Shell, 1369f4259b30SLisandro Dalcin /*114*/ NULL, 1370f4259b30SLisandro Dalcin NULL, 1371f4259b30SLisandro Dalcin NULL, 1372f4259b30SLisandro Dalcin NULL, 1373f4259b30SLisandro Dalcin NULL, 1374f4259b30SLisandro Dalcin /*119*/ NULL, 1375f4259b30SLisandro Dalcin NULL, 1376f4259b30SLisandro Dalcin NULL, 1377f4259b30SLisandro Dalcin NULL, 1378f4259b30SLisandro Dalcin NULL, 1379f4259b30SLisandro Dalcin /*124*/ NULL, 1380f4259b30SLisandro Dalcin NULL, 1381f4259b30SLisandro Dalcin NULL, 1382f4259b30SLisandro Dalcin NULL, 1383f4259b30SLisandro Dalcin NULL, 1384f4259b30SLisandro Dalcin /*129*/ NULL, 1385f4259b30SLisandro Dalcin NULL, 1386f4259b30SLisandro Dalcin NULL, 1387f4259b30SLisandro Dalcin NULL, 1388f4259b30SLisandro Dalcin NULL, 1389f4259b30SLisandro Dalcin /*134*/ NULL, 1390f4259b30SLisandro Dalcin NULL, 1391f4259b30SLisandro Dalcin NULL, 1392f4259b30SLisandro Dalcin NULL, 1393f4259b30SLisandro Dalcin NULL, 1394f4259b30SLisandro Dalcin /*139*/ NULL, 1395f4259b30SLisandro Dalcin NULL, 1396d70f29a3SPierre Jolivet NULL, 1397d70f29a3SPierre Jolivet NULL, 1398d70f29a3SPierre Jolivet NULL, 1399d70f29a3SPierre Jolivet /*144*/ NULL, 1400d70f29a3SPierre Jolivet NULL, 1401d70f29a3SPierre Jolivet NULL, 140299a7f59eSMark Adams NULL, 140399a7f59eSMark Adams NULL, 14047fb60732SBarry Smith NULL, 14059371c9d4SSatish Balay /*150*/ NULL}; 1406273d9f13SBarry Smith 14079371c9d4SSatish Balay static PetscErrorCode MatShellSetContext_Shell(Mat mat, void *ctx) { 1408789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell *)mat->data; 1409789d8953SBarry Smith 1410789d8953SBarry Smith PetscFunctionBegin; 1411800f99ffSJeremy L Thompson if (ctx) { 1412800f99ffSJeremy L Thompson PetscContainer ctxcontainer; 1413800f99ffSJeremy L Thompson PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)mat), &ctxcontainer)); 1414800f99ffSJeremy L Thompson PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 1415800f99ffSJeremy L Thompson PetscCall(PetscObjectCompose((PetscObject)mat, "MatShell ctx", (PetscObject)ctxcontainer)); 1416800f99ffSJeremy L Thompson shell->ctxcontainer = ctxcontainer; 1417800f99ffSJeremy L Thompson PetscCall(PetscContainerDestroy(&ctxcontainer)); 1418800f99ffSJeremy L Thompson } else { 1419800f99ffSJeremy L Thompson PetscCall(PetscObjectCompose((PetscObject)mat, "MatShell ctx", NULL)); 1420800f99ffSJeremy L Thompson shell->ctxcontainer = NULL; 1421800f99ffSJeremy L Thompson } 1422800f99ffSJeremy L Thompson 1423800f99ffSJeremy L Thompson PetscFunctionReturn(0); 1424800f99ffSJeremy L Thompson } 1425800f99ffSJeremy L Thompson 14269371c9d4SSatish Balay PetscErrorCode MatShellSetContextDestroy_Shell(Mat mat, PetscErrorCode (*f)(void *)) { 1427800f99ffSJeremy L Thompson Mat_Shell *shell = (Mat_Shell *)mat->data; 1428800f99ffSJeremy L Thompson 1429800f99ffSJeremy L Thompson PetscFunctionBegin; 1430800f99ffSJeremy L Thompson if (shell->ctxcontainer) PetscCall(PetscContainerSetUserDestroy(shell->ctxcontainer, f)); 1431789d8953SBarry Smith PetscFunctionReturn(0); 1432789d8953SBarry Smith } 1433789d8953SBarry Smith 14349371c9d4SSatish Balay static PetscErrorCode MatShellSetVecType_Shell(Mat mat, VecType vtype) { 1435db77db73SJeremy L Thompson PetscFunctionBegin; 14369566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->defaultvectype)); 14379566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(vtype, (char **)&mat->defaultvectype)); 1438db77db73SJeremy L Thompson PetscFunctionReturn(0); 1439db77db73SJeremy L Thompson } 1440db77db73SJeremy L Thompson 14419371c9d4SSatish Balay PetscErrorCode MatShellSetManageScalingShifts_Shell(Mat A) { 1442789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell *)A->data; 1443789d8953SBarry Smith 1444789d8953SBarry Smith PetscFunctionBegin; 1445789d8953SBarry Smith shell->managescalingshifts = PETSC_FALSE; 1446789d8953SBarry Smith A->ops->diagonalset = NULL; 1447789d8953SBarry Smith A->ops->diagonalscale = NULL; 1448789d8953SBarry Smith A->ops->scale = NULL; 1449789d8953SBarry Smith A->ops->shift = NULL; 1450789d8953SBarry Smith A->ops->axpy = NULL; 1451789d8953SBarry Smith PetscFunctionReturn(0); 1452789d8953SBarry Smith } 1453789d8953SBarry Smith 14549371c9d4SSatish Balay static PetscErrorCode MatShellSetOperation_Shell(Mat mat, MatOperation op, void (*f)(void)) { 1455feb237baSPierre Jolivet Mat_Shell *shell = (Mat_Shell *)mat->data; 1456789d8953SBarry Smith 1457789d8953SBarry Smith PetscFunctionBegin; 1458789d8953SBarry Smith switch (op) { 14599371c9d4SSatish Balay case MATOP_DESTROY: shell->ops->destroy = (PetscErrorCode(*)(Mat))f; break; 1460789d8953SBarry Smith case MATOP_VIEW: 1461ad540459SPierre Jolivet if (!mat->ops->viewnative) mat->ops->viewnative = mat->ops->view; 1462789d8953SBarry Smith mat->ops->view = (PetscErrorCode(*)(Mat, PetscViewer))f; 1463789d8953SBarry Smith break; 14649371c9d4SSatish Balay case MATOP_COPY: shell->ops->copy = (PetscErrorCode(*)(Mat, Mat, MatStructure))f; break; 1465789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1466789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1467789d8953SBarry Smith case MATOP_SHIFT: 1468789d8953SBarry Smith case MATOP_SCALE: 1469789d8953SBarry Smith case MATOP_AXPY: 1470789d8953SBarry Smith case MATOP_ZERO_ROWS: 1471789d8953SBarry Smith case MATOP_ZERO_ROWS_COLUMNS: 147228b400f6SJacob Faibussowitsch PetscCheck(!shell->managescalingshifts, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()"); 1473789d8953SBarry Smith (((void (**)(void))mat->ops)[op]) = f; 1474789d8953SBarry Smith break; 1475789d8953SBarry Smith case MATOP_GET_DIAGONAL: 1476789d8953SBarry Smith if (shell->managescalingshifts) { 1477789d8953SBarry Smith shell->ops->getdiagonal = (PetscErrorCode(*)(Mat, Vec))f; 1478789d8953SBarry Smith mat->ops->getdiagonal = MatGetDiagonal_Shell; 1479789d8953SBarry Smith } else { 1480789d8953SBarry Smith shell->ops->getdiagonal = NULL; 1481789d8953SBarry Smith mat->ops->getdiagonal = (PetscErrorCode(*)(Mat, Vec))f; 1482789d8953SBarry Smith } 1483789d8953SBarry Smith break; 1484789d8953SBarry Smith case MATOP_MULT: 1485789d8953SBarry Smith if (shell->managescalingshifts) { 1486789d8953SBarry Smith shell->ops->mult = (PetscErrorCode(*)(Mat, Vec, Vec))f; 1487789d8953SBarry Smith mat->ops->mult = MatMult_Shell; 1488789d8953SBarry Smith } else { 1489789d8953SBarry Smith shell->ops->mult = NULL; 1490789d8953SBarry Smith mat->ops->mult = (PetscErrorCode(*)(Mat, Vec, Vec))f; 1491789d8953SBarry Smith } 1492789d8953SBarry Smith break; 1493789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 1494789d8953SBarry Smith if (shell->managescalingshifts) { 1495789d8953SBarry Smith shell->ops->multtranspose = (PetscErrorCode(*)(Mat, Vec, Vec))f; 1496789d8953SBarry Smith mat->ops->multtranspose = MatMultTranspose_Shell; 1497789d8953SBarry Smith } else { 1498789d8953SBarry Smith shell->ops->multtranspose = NULL; 1499789d8953SBarry Smith mat->ops->multtranspose = (PetscErrorCode(*)(Mat, Vec, Vec))f; 1500789d8953SBarry Smith } 1501789d8953SBarry Smith break; 15029371c9d4SSatish Balay default: (((void (**)(void))mat->ops)[op]) = f; break; 1503789d8953SBarry Smith } 1504789d8953SBarry Smith PetscFunctionReturn(0); 1505789d8953SBarry Smith } 1506789d8953SBarry Smith 15079371c9d4SSatish Balay static PetscErrorCode MatShellGetOperation_Shell(Mat mat, MatOperation op, void (**f)(void)) { 1508789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell *)mat->data; 1509789d8953SBarry Smith 1510789d8953SBarry Smith PetscFunctionBegin; 1511789d8953SBarry Smith switch (op) { 15129371c9d4SSatish Balay case MATOP_DESTROY: *f = (void (*)(void))shell->ops->destroy; break; 15139371c9d4SSatish Balay case MATOP_VIEW: *f = (void (*)(void))mat->ops->view; break; 15149371c9d4SSatish Balay case MATOP_COPY: *f = (void (*)(void))shell->ops->copy; break; 1515789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1516789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1517789d8953SBarry Smith case MATOP_SHIFT: 1518789d8953SBarry Smith case MATOP_SCALE: 1519789d8953SBarry Smith case MATOP_AXPY: 1520789d8953SBarry Smith case MATOP_ZERO_ROWS: 15219371c9d4SSatish Balay case MATOP_ZERO_ROWS_COLUMNS: *f = (((void (**)(void))mat->ops)[op]); break; 1522789d8953SBarry Smith case MATOP_GET_DIAGONAL: 15239371c9d4SSatish Balay if (shell->ops->getdiagonal) *f = (void (*)(void))shell->ops->getdiagonal; 15249371c9d4SSatish Balay else *f = (((void (**)(void))mat->ops)[op]); 1525789d8953SBarry Smith break; 1526789d8953SBarry Smith case MATOP_MULT: 15279371c9d4SSatish Balay if (shell->ops->mult) *f = (void (*)(void))shell->ops->mult; 15289371c9d4SSatish Balay else *f = (((void (**)(void))mat->ops)[op]); 1529789d8953SBarry Smith break; 1530789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 15319371c9d4SSatish Balay if (shell->ops->multtranspose) *f = (void (*)(void))shell->ops->multtranspose; 15329371c9d4SSatish Balay else *f = (((void (**)(void))mat->ops)[op]); 1533789d8953SBarry Smith break; 15349371c9d4SSatish Balay default: *f = (((void (**)(void))mat->ops)[op]); 1535789d8953SBarry Smith } 1536789d8953SBarry Smith PetscFunctionReturn(0); 1537789d8953SBarry Smith } 1538789d8953SBarry Smith 15390bad9183SKris Buschelman /*MC 1540fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 15410bad9183SKris Buschelman 15420bad9183SKris Buschelman Level: advanced 15430bad9183SKris Buschelman 154411a5261eSBarry Smith .seealso: `Mat`, `MatCreateShell()` 15450bad9183SKris Buschelman M*/ 15460bad9183SKris Buschelman 15479371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) { 1548273d9f13SBarry Smith Mat_Shell *b; 1549273d9f13SBarry Smith 1550273d9f13SBarry Smith PetscFunctionBegin; 15519566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A->ops, &MatOps_Values, sizeof(struct _MatOps))); 1552273d9f13SBarry Smith 1553*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 1554273d9f13SBarry Smith A->data = (void *)b; 1555273d9f13SBarry Smith 1556800f99ffSJeremy L Thompson b->ctxcontainer = NULL; 1557ef66eb69SBarry Smith b->vshift = 0.0; 1558ef66eb69SBarry Smith b->vscale = 1.0; 15590c0fd78eSBarry Smith b->managescalingshifts = PETSC_TRUE; 1560273d9f13SBarry Smith A->assembled = PETSC_TRUE; 1561210f0121SBarry Smith A->preallocated = PETSC_FALSE; 15622205254eSKarl Rupp 15639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellGetContext_C", MatShellGetContext_Shell)); 15649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Shell)); 1565800f99ffSJeremy L Thompson PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Shell)); 15669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetVecType_C", MatShellSetVecType_Shell)); 15679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Shell)); 15689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetOperation_C", MatShellSetOperation_Shell)); 15699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellGetOperation_C", MatShellGetOperation_Shell)); 15709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetMatProductOperation_C", MatShellSetMatProductOperation_Shell)); 15719566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATSHELL)); 1572273d9f13SBarry Smith PetscFunctionReturn(0); 1573273d9f13SBarry Smith } 1574e51e0e81SBarry Smith 15754b828684SBarry Smith /*@C 157611a5261eSBarry Smith MatCreateShell - Creates a new matrix of `MatType` `MATSHELL` for use with a user-defined 1577ff756334SLois Curfman McInnes private data storage format. 1578e51e0e81SBarry Smith 1579d083f849SBarry Smith Collective 1580c7fcc2eaSBarry Smith 1581e51e0e81SBarry Smith Input Parameters: 1582c7fcc2eaSBarry Smith + comm - MPI communicator 1583c7fcc2eaSBarry Smith . m - number of local rows (must be given) 1584c7fcc2eaSBarry Smith . n - number of local columns (must be given) 1585c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 1586c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 1587c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 1588e51e0e81SBarry Smith 1589ff756334SLois Curfman McInnes Output Parameter: 159044cd7ae7SLois Curfman McInnes . A - the matrix 1591e51e0e81SBarry Smith 1592ff2fd236SBarry Smith Level: advanced 1593ff2fd236SBarry Smith 1594f39d1f56SLois Curfman McInnes Usage: 15955bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1596f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 1597c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1598f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 1599f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 1600f39d1f56SLois Curfman McInnes 1601ff756334SLois Curfman McInnes Notes: 1602ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 160311a5261eSBarry Smith with `KSP` (such as, for use with matrix-free methods). You should not 1604ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 1605e51e0e81SBarry Smith 1606f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 1607f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 1608645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 160911a5261eSBarry Smith products using `MatMult`(A,x,y) the user should set the number 1610645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 1611645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 1612645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 1613645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 1614645985a0SLois Curfman McInnes For example, 1615f39d1f56SLois Curfman McInnes 1616f39d1f56SLois Curfman McInnes $ 1617f39d1f56SLois Curfman McInnes $ Vec x, y 16185bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1619645985a0SLois Curfman McInnes $ Mat A 1620f39d1f56SLois Curfman McInnes $ 1621c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 1622c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 1623f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 1624c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 1625c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 1626c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1627645985a0SLois Curfman McInnes $ MatMult(A,x,y); 162845960306SStefano Zampini $ MatDestroy(&A); 162945960306SStefano Zampini $ VecDestroy(&y); 163045960306SStefano Zampini $ VecDestroy(&x); 1631645985a0SLois Curfman McInnes $ 1632e51e0e81SBarry Smith 163311a5261eSBarry Smith `MATSHELL` handles `MatShift()`, `MatDiagonalSet()`, `MatDiagonalScale()`, `MatAXPY()`, `MatScale()`, `MatZeroRows()` and `MatZeroRowsColumns()` internally, so these 163411a5261eSBarry Smith operations cannot be overwritten unless `MatShellSetManageScalingShifts()` is called. 16359b53d762SBarry Smith 16360c0fd78eSBarry Smith For rectangular matrices do all the scalings and shifts make sense? 16370c0fd78eSBarry Smith 163895452b02SPatrick Sanan Developers Notes: 163995452b02SPatrick Sanan Regarding shifting and scaling. The general form is 16400c0fd78eSBarry Smith 16410c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 16420c0fd78eSBarry Smith 16430c0fd78eSBarry Smith The order you apply the operations is important. For example if you have a dshift then 16440c0fd78eSBarry Smith apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift 16450c0fd78eSBarry Smith you get s*vscale*A + diag(shift) 16460c0fd78eSBarry Smith 16470c0fd78eSBarry Smith A is the user provided function. 16480c0fd78eSBarry Smith 164911a5261eSBarry Smith `KSP`/`PC` uses changes in the` Mat`'s "state" to decide if preconditioners need to be rebuilt: `PCSetUp()` only calls the setup() for 165011a5261eSBarry Smith for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a `MATSHELL` to trigger 165111a5261eSBarry Smith an update in the preconditioner you must call `MatAssemblyBegin()` and `MatAssemblyEnd()` or `PetscObjectStateIncrease`((`PetscObject`)mat); 165211a5261eSBarry Smith each time the `MATSHELL` matrix has changed. 1653ad2f5c3fSBarry Smith 165411a5261eSBarry Smith Matrix product operations (i.e. `MatMat()`, `MatTransposeMat()` etc) can be specified using `MatShellSetMatProductOperation()` 1655b77ba244SStefano Zampini 165611a5261eSBarry Smith Calling `MatAssemblyBegin()`/`MatAssemblyEnd()` on a `MATSHELL` removes any previously supplied shift and scales that were provided 165711a5261eSBarry Smith with `MatDiagonalSet()`, `MatShift()`, `MatScale()`, or `MatDiagonalScale()`. 1658ad2f5c3fSBarry Smith 165911a5261eSBarry Smith Fortran Note: 166011a5261eSBarry Smith To use this from Fortran with a ctx you must write an interface definition for this 166111a5261eSBarry Smith function and for `MatShellGetContext()` that tells Fortran the Fortran derived data type you are passing 166211a5261eSBarry Smith in as the ctx argument. 166311a5261eSBarry Smith 166411a5261eSBarry Smith .seealso: `MATSHELL`, `MatShellSetOperation()`, `MatHasOperation()`, `MatShellGetContext()`, `MatShellSetContext()`, `MATSHELL`, `MatShellSetManageScalingShifts()`, `MatShellSetMatProductOperation()` 1665e51e0e81SBarry Smith @*/ 16669371c9d4SSatish Balay PetscErrorCode MatCreateShell(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, void *ctx, Mat *A) { 16673a40ed3dSBarry Smith PetscFunctionBegin; 16689566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, A)); 16699566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A, m, n, M, N)); 16709566063dSJacob Faibussowitsch PetscCall(MatSetType(*A, MATSHELL)); 16719566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(*A, ctx)); 16729566063dSJacob Faibussowitsch PetscCall(MatSetUp(*A)); 1673273d9f13SBarry Smith PetscFunctionReturn(0); 1674c7fcc2eaSBarry Smith } 1675c7fcc2eaSBarry Smith 1676c6866cfdSSatish Balay /*@ 167711a5261eSBarry Smith MatShellSetContext - sets the context for a `MATSHELL` shell matrix 1678c7fcc2eaSBarry Smith 167911a5261eSBarry Smith Logically Collective on mat 1680c7fcc2eaSBarry Smith 1681273d9f13SBarry Smith Input Parameters: 168211a5261eSBarry Smith + mat - the `MATSHELL` shell matrix 1683273d9f13SBarry Smith - ctx - the context 1684273d9f13SBarry Smith 1685273d9f13SBarry Smith Level: advanced 1686273d9f13SBarry Smith 168711a5261eSBarry Smith Fortran Note: 168895452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 1689daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 1690273d9f13SBarry Smith 169111a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()` 16920bc0a719SBarry Smith @*/ 16939371c9d4SSatish Balay PetscErrorCode MatShellSetContext(Mat mat, void *ctx) { 1694273d9f13SBarry Smith PetscFunctionBegin; 16950700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1696cac4c232SBarry Smith PetscTryMethod(mat, "MatShellSetContext_C", (Mat, void *), (mat, ctx)); 16973a40ed3dSBarry Smith PetscFunctionReturn(0); 1698e51e0e81SBarry Smith } 1699e51e0e81SBarry Smith 1700800f99ffSJeremy L Thompson /*@ 170111a5261eSBarry Smith MatShellSetContextDestroy - sets the destroy function for a `MATSHELL` shell matrix context 1702800f99ffSJeremy L Thompson 170311a5261eSBarry Smith Logically Collective on mat 1704800f99ffSJeremy L Thompson 1705800f99ffSJeremy L Thompson Input Parameters: 1706800f99ffSJeremy L Thompson + mat - the shell matrix 1707800f99ffSJeremy L Thompson - f - the context destroy function 1708800f99ffSJeremy L Thompson 1709800f99ffSJeremy L Thompson Note: 1710800f99ffSJeremy L Thompson If the `MatShell` is never duplicated, the behavior of this function is equivalent 171111a5261eSBarry Smith to `MatShellSetOperation`(`Mat`,`MATOP_DESTROY`,f). However, `MatShellSetContextDestroy()` 1712800f99ffSJeremy L Thompson ensures proper reference counting for the user provided context data in the case that 171311a5261eSBarry Smith the `MATSHELL` is duplicated. 1714800f99ffSJeremy L Thompson 1715800f99ffSJeremy L Thompson Level: advanced 1716800f99ffSJeremy L Thompson 171711a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateShell()`, `MatShellSetContext()` 1718800f99ffSJeremy L Thompson @*/ 17199371c9d4SSatish Balay PetscErrorCode MatShellSetContextDestroy(Mat mat, PetscErrorCode (*f)(void *)) { 1720800f99ffSJeremy L Thompson PetscFunctionBegin; 1721800f99ffSJeremy L Thompson PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1722800f99ffSJeremy L Thompson PetscTryMethod(mat, "MatShellSetContextDestroy_C", (Mat, PetscErrorCode(*)(void *)), (mat, f)); 1723800f99ffSJeremy L Thompson PetscFunctionReturn(0); 1724800f99ffSJeremy L Thompson } 1725800f99ffSJeremy L Thompson 1726db77db73SJeremy L Thompson /*@C 172711a5261eSBarry Smith MatShellSetVecType - Sets the type of `Vec` returned by `MatCreateVecs()` 1728db77db73SJeremy L Thompson 1729db77db73SJeremy L Thompson Logically collective 1730db77db73SJeremy L Thompson 1731db77db73SJeremy L Thompson Input Parameters: 173211a5261eSBarry Smith + mat - the `MATSHELL` shell matrix 1733db77db73SJeremy L Thompson - vtype - type to use for creating vectors 1734db77db73SJeremy L Thompson 1735db77db73SJeremy L Thompson Level: advanced 1736db77db73SJeremy L Thompson 173711a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateVecs()` 1738db77db73SJeremy L Thompson @*/ 17399371c9d4SSatish Balay PetscErrorCode MatShellSetVecType(Mat mat, VecType vtype) { 1740db77db73SJeremy L Thompson PetscFunctionBegin; 1741cac4c232SBarry Smith PetscTryMethod(mat, "MatShellSetVecType_C", (Mat, VecType), (mat, vtype)); 1742db77db73SJeremy L Thompson PetscFunctionReturn(0); 1743db77db73SJeremy L Thompson } 1744db77db73SJeremy L Thompson 17450c0fd78eSBarry Smith /*@ 174611a5261eSBarry Smith MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the `MATSHELL`. Must be called immediately 174711a5261eSBarry Smith after `MatCreateShell()` 17480c0fd78eSBarry Smith 174911a5261eSBarry Smith Logically Collective on A 17500c0fd78eSBarry Smith 17510c0fd78eSBarry Smith Input Parameter: 175211a5261eSBarry Smith . mat - the `MATSHELL` shell matrix 17530c0fd78eSBarry Smith 17540c0fd78eSBarry Smith Level: advanced 17550c0fd78eSBarry Smith 175611a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellSetContext()`, `MatShellSetOperation()` 17570c0fd78eSBarry Smith @*/ 17589371c9d4SSatish Balay PetscErrorCode MatShellSetManageScalingShifts(Mat A) { 17590c0fd78eSBarry Smith PetscFunctionBegin; 17600c0fd78eSBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 1761cac4c232SBarry Smith PetscTryMethod(A, "MatShellSetManageScalingShifts_C", (Mat), (A)); 17620c0fd78eSBarry Smith PetscFunctionReturn(0); 17630c0fd78eSBarry Smith } 17640c0fd78eSBarry Smith 1765c16cb8f2SBarry Smith /*@C 176611a5261eSBarry Smith MatShellTestMult - Compares the multiply routine provided to the `MATSHELL` with differencing on a given function. 1767f3b1f45cSBarry Smith 176811a5261eSBarry Smith Logically Collective on mat 1769f3b1f45cSBarry Smith 1770f3b1f45cSBarry Smith Input Parameters: 177111a5261eSBarry Smith + mat - the `MATSHELL` shell matrix 1772f3b1f45cSBarry Smith . f - the function 177311a5261eSBarry Smith . base - differences are computed around this vector, see `MatMFFDSetBase()`, for Jacobians this is the point at which the Jacobian is being evaluated 1774f3b1f45cSBarry Smith - ctx - an optional context for the function 1775f3b1f45cSBarry Smith 1776f3b1f45cSBarry Smith Output Parameter: 177711a5261eSBarry Smith . flg - `PETSC_TRUE` if the multiply is likely correct 1778f3b1f45cSBarry Smith 1779f3b1f45cSBarry Smith Options Database: 1780f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1781f3b1f45cSBarry Smith 1782f3b1f45cSBarry Smith Level: advanced 1783f3b1f45cSBarry Smith 178411a5261eSBarry Smith Fortran Note: 178595452b02SPatrick Sanan Not supported from Fortran 1786f3b1f45cSBarry Smith 178711a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellTestMultTranspose()` 1788f3b1f45cSBarry Smith @*/ 17899371c9d4SSatish Balay PetscErrorCode MatShellTestMult(Mat mat, PetscErrorCode (*f)(void *, Vec, Vec), Vec base, void *ctx, PetscBool *flg) { 1790f3b1f45cSBarry Smith PetscInt m, n; 1791f3b1f45cSBarry Smith Mat mf, Dmf, Dmat, Ddiff; 1792f3b1f45cSBarry Smith PetscReal Diffnorm, Dmfnorm; 179374e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1794f3b1f45cSBarry Smith 1795f3b1f45cSBarry Smith PetscFunctionBegin; 1796f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 17979566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, ((PetscObject)mat)->prefix, "-mat_shell_test_mult_view", &v)); 17989566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat, &m, &n)); 17999566063dSJacob Faibussowitsch PetscCall(MatCreateMFFD(PetscObjectComm((PetscObject)mat), m, n, PETSC_DECIDE, PETSC_DECIDE, &mf)); 18009566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(mf, f, ctx)); 18019566063dSJacob Faibussowitsch PetscCall(MatMFFDSetBase(mf, base, NULL)); 1802f3b1f45cSBarry Smith 18039566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(mf, MATAIJ, &Dmf)); 18049566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(mat, MATAIJ, &Dmat)); 1805f3b1f45cSBarry Smith 18069566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Dmat, MAT_COPY_VALUES, &Ddiff)); 18079566063dSJacob Faibussowitsch PetscCall(MatAXPY(Ddiff, -1.0, Dmf, DIFFERENT_NONZERO_PATTERN)); 18089566063dSJacob Faibussowitsch PetscCall(MatNorm(Ddiff, NORM_FROBENIUS, &Diffnorm)); 18099566063dSJacob Faibussowitsch PetscCall(MatNorm(Dmf, NORM_FROBENIUS, &Dmfnorm)); 1810f3b1f45cSBarry Smith if (Diffnorm / Dmfnorm > 10 * PETSC_SQRT_MACHINE_EPSILON) { 1811f3b1f45cSBarry Smith flag = PETSC_FALSE; 1812f3b1f45cSBarry Smith if (v) { 18139566063dSJacob Faibussowitsch PetscCall(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))); 18149566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Ddiff, (PetscObject)mat, "-mat_shell_test_mult_view")); 18159566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Dmf, (PetscObject)mat, "-mat_shell_test_mult_view")); 18169566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Dmat, (PetscObject)mat, "-mat_shell_test_mult_view")); 1817f3b1f45cSBarry Smith } 1818f3b1f45cSBarry Smith } else if (v) { 18199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat), "MATSHELL and matrix free multiple appear to produce the same results\n")); 1820f3b1f45cSBarry Smith } 1821f3b1f45cSBarry Smith if (flg) *flg = flag; 18229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ddiff)); 18239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mf)); 18249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Dmf)); 18259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Dmat)); 1826f3b1f45cSBarry Smith PetscFunctionReturn(0); 1827f3b1f45cSBarry Smith } 1828f3b1f45cSBarry Smith 1829f3b1f45cSBarry Smith /*@C 183011a5261eSBarry Smith MatShellTestMultTranspose - Compares the multiply transpose routine provided to the `MATSHELL` with differencing on a given function. 1831f3b1f45cSBarry Smith 183211a5261eSBarry Smith Logically Collective on mat 1833f3b1f45cSBarry Smith 1834f3b1f45cSBarry Smith Input Parameters: 183511a5261eSBarry Smith + mat - the `MATSHELL` shell matrix 1836f3b1f45cSBarry Smith . f - the function 183711a5261eSBarry Smith . base - differences are computed around this vector, see `MatMFFDSetBase()`, for Jacobians this is the point at which the Jacobian is being evaluated 1838f3b1f45cSBarry Smith - ctx - an optional context for the function 1839f3b1f45cSBarry Smith 1840f3b1f45cSBarry Smith Output Parameter: 184111a5261eSBarry Smith . flg - `PETSC_TRUE` if the multiply is likely correct 1842f3b1f45cSBarry Smith 1843f3b1f45cSBarry Smith Options Database: 1844f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1845f3b1f45cSBarry Smith 1846f3b1f45cSBarry Smith Level: advanced 1847f3b1f45cSBarry Smith 184811a5261eSBarry Smith Fortran Note: 184995452b02SPatrick Sanan Not supported from Fortran 1850f3b1f45cSBarry Smith 185111a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellTestMult()` 1852f3b1f45cSBarry Smith @*/ 18539371c9d4SSatish Balay PetscErrorCode MatShellTestMultTranspose(Mat mat, PetscErrorCode (*f)(void *, Vec, Vec), Vec base, void *ctx, PetscBool *flg) { 1854f3b1f45cSBarry Smith Vec x, y, z; 1855f3b1f45cSBarry Smith PetscInt m, n, M, N; 1856f3b1f45cSBarry Smith Mat mf, Dmf, Dmat, Ddiff; 1857f3b1f45cSBarry Smith PetscReal Diffnorm, Dmfnorm; 185874e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1859f3b1f45cSBarry Smith 1860f3b1f45cSBarry Smith PetscFunctionBegin; 1861f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 18629566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, ((PetscObject)mat)->prefix, "-mat_shell_test_mult_transpose_view", &v)); 18639566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat, &x, &y)); 18649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &z)); 18659566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat, &m, &n)); 18669566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat, &M, &N)); 18679566063dSJacob Faibussowitsch PetscCall(MatCreateMFFD(PetscObjectComm((PetscObject)mat), m, n, M, N, &mf)); 18689566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(mf, f, ctx)); 18699566063dSJacob Faibussowitsch PetscCall(MatMFFDSetBase(mf, base, NULL)); 18709566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(mf, MATAIJ, &Dmf)); 18719566063dSJacob Faibussowitsch PetscCall(MatTranspose(Dmf, MAT_INPLACE_MATRIX, &Dmf)); 18729566063dSJacob Faibussowitsch PetscCall(MatComputeOperatorTranspose(mat, MATAIJ, &Dmat)); 1873f3b1f45cSBarry Smith 18749566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Dmat, MAT_COPY_VALUES, &Ddiff)); 18759566063dSJacob Faibussowitsch PetscCall(MatAXPY(Ddiff, -1.0, Dmf, DIFFERENT_NONZERO_PATTERN)); 18769566063dSJacob Faibussowitsch PetscCall(MatNorm(Ddiff, NORM_FROBENIUS, &Diffnorm)); 18779566063dSJacob Faibussowitsch PetscCall(MatNorm(Dmf, NORM_FROBENIUS, &Dmfnorm)); 1878f3b1f45cSBarry Smith if (Diffnorm / Dmfnorm > 10 * PETSC_SQRT_MACHINE_EPSILON) { 1879f3b1f45cSBarry Smith flag = PETSC_FALSE; 1880f3b1f45cSBarry Smith if (v) { 18819566063dSJacob Faibussowitsch PetscCall(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))); 18829566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Ddiff, (PetscObject)mat, "-mat_shell_test_mult_transpose_view")); 18839566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Dmf, (PetscObject)mat, "-mat_shell_test_mult_transpose_view")); 18849566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Dmat, (PetscObject)mat, "-mat_shell_test_mult_transpose_view")); 1885f3b1f45cSBarry Smith } 1886f3b1f45cSBarry Smith } else if (v) { 18879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat), "MATSHELL transpose and matrix free multiple appear to produce the same results\n")); 1888f3b1f45cSBarry Smith } 1889f3b1f45cSBarry Smith if (flg) *flg = flag; 18909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mf)); 18919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Dmat)); 18929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ddiff)); 18939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Dmf)); 18949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 18959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 18969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 1897f3b1f45cSBarry Smith PetscFunctionReturn(0); 1898f3b1f45cSBarry Smith } 1899f3b1f45cSBarry Smith 1900f3b1f45cSBarry Smith /*@C 190111a5261eSBarry Smith MatShellSetOperation - Allows user to set a matrix operation for a `MATSHELL` shell matrix. 1902e51e0e81SBarry Smith 190311a5261eSBarry Smith Logically Collective on mat 1904fee21e36SBarry Smith 1905c7fcc2eaSBarry Smith Input Parameters: 190611a5261eSBarry Smith + mat - the `MATSHELL` shell matrix 1907c7fcc2eaSBarry Smith . op - the name of the operation 1908789d8953SBarry Smith - g - the function that provides the operation. 1909c7fcc2eaSBarry Smith 191015091d37SBarry Smith Level: advanced 191115091d37SBarry Smith 1912fae171e0SBarry Smith Usage: 1913e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 1914b77ba244SStefano Zampini $ MatCreateShell(comm,m,n,M,N,ctx,&A); 1915b77ba244SStefano Zampini $ MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 19160b627109SLois Curfman McInnes 1917a62d957aSLois Curfman McInnes Notes: 1918e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 19191c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 1920a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 192111a5261eSBarry Smith user interface routine (e.g., `MatMult()` -> `MATOP_MULT`). 1922a62d957aSLois Curfman McInnes 192311a5261eSBarry Smith All user-provided functions (except for `MATOP_DESTROY`) should have the same calling 1924deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 1925deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 1926deebb3c3SLois Curfman McInnes routines, e.g., 1927a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1928a62d957aSLois Curfman McInnes 19294aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 19304aa34b0aSBarry Smith nonzero on failure. 19314aa34b0aSBarry Smith 1932a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 193311a5261eSBarry Smith `MatShellGetContext()` to obtain the user-defined context that was 193411a5261eSBarry Smith set by `MatCreateShell()`. 1935a62d957aSLois Curfman McInnes 193611a5261eSBarry Smith Use `MatSetOperation()` to set an operation for any matrix type. For matrix product operations (i.e. `MatMatXXX()`, `MatTransposeMatXXX()` etc) use `MatShellSetMatProductOperation()` 1937b77ba244SStefano Zampini 193811a5261eSBarry Smith Fortran Note: 193911a5261eSBarry Smith For `MatCreateVecs()` the user code should check if the input left or right matrix is -1 and in that case not 1940c4762a1bSJed Brown generate a matrix. See src/mat/tests/ex120f.F 1941501d9185SBarry Smith 194211a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellSetContext()`, `MatSetOperation()`, `MatShellSetManageScalingShifts()`, `MatShellSetMatProductOperation()` 1943e51e0e81SBarry Smith @*/ 19449371c9d4SSatish Balay PetscErrorCode MatShellSetOperation(Mat mat, MatOperation op, void (*g)(void)) { 19453a40ed3dSBarry Smith PetscFunctionBegin; 19460700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1947cac4c232SBarry Smith PetscTryMethod(mat, "MatShellSetOperation_C", (Mat, MatOperation, void (*)(void)), (mat, op, g)); 19483a40ed3dSBarry Smith PetscFunctionReturn(0); 1949e51e0e81SBarry Smith } 1950f0479e8cSBarry Smith 1951d4bb536fSBarry Smith /*@C 195211a5261eSBarry Smith MatShellGetOperation - Gets a matrix function for a `MATSHELL` shell matrix. 1953d4bb536fSBarry Smith 1954c7fcc2eaSBarry Smith Not Collective 1955c7fcc2eaSBarry Smith 1956d4bb536fSBarry Smith Input Parameters: 195711a5261eSBarry Smith + mat - the `MATSHELL` shell matrix 1958c7fcc2eaSBarry Smith - op - the name of the operation 1959d4bb536fSBarry Smith 1960d4bb536fSBarry Smith Output Parameter: 1961789d8953SBarry Smith . g - the function that provides the operation. 1962d4bb536fSBarry Smith 196315091d37SBarry Smith Level: advanced 196415091d37SBarry Smith 196511a5261eSBarry Smith Note: 1966e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 1967d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 1968d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 196911a5261eSBarry Smith user interface routine (e.g., `MatMult()` -> `MATOP_MULT`). 1970d4bb536fSBarry Smith 1971d4bb536fSBarry Smith All user-provided functions have the same calling 1972d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 1973d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 1974d4bb536fSBarry Smith routines, e.g., 1975d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1976d4bb536fSBarry Smith 1977d4bb536fSBarry Smith Within each user-defined routine, the user should call 197811a5261eSBarry Smith `MatShellGetContext()` to obtain the user-defined context that was 197911a5261eSBarry Smith set by `MatCreateShell()`. 1980d4bb536fSBarry Smith 198111a5261eSBarry Smith .seealso: `MATSHELL`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellSetOperation()`, `MatShellSetContext()` 1982d4bb536fSBarry Smith @*/ 19839371c9d4SSatish Balay PetscErrorCode MatShellGetOperation(Mat mat, MatOperation op, void (**g)(void)) { 19843a40ed3dSBarry Smith PetscFunctionBegin; 19850700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1986cac4c232SBarry Smith PetscUseMethod(mat, "MatShellGetOperation_C", (Mat, MatOperation, void (**)(void)), (mat, op, g)); 19873a40ed3dSBarry Smith PetscFunctionReturn(0); 1988d4bb536fSBarry Smith } 1989b77ba244SStefano Zampini 1990b77ba244SStefano Zampini /*@ 199111a5261eSBarry Smith MatIsShell - Inquires if a matrix is derived from `MATSHELL` 1992b77ba244SStefano Zampini 1993b77ba244SStefano Zampini Input Parameter: 1994b77ba244SStefano Zampini . mat - the matrix 1995b77ba244SStefano Zampini 1996b77ba244SStefano Zampini Output Parameter: 1997b77ba244SStefano Zampini . flg - the boolean value 1998b77ba244SStefano Zampini 1999b77ba244SStefano Zampini Level: developer 2000b77ba244SStefano Zampini 200111a5261eSBarry Smith Developer Note: 2002013e2dc7SBarry Smith In the future, we should allow the object type name to be changed still using the `MATSHELL` data structure for other matrices (i.e. `MATTRANSPOSEVIRTUAL`, `MATSCHURCOMPLEMENT` etc) 2003b77ba244SStefano Zampini 2004013e2dc7SBarry Smith .seealso: `MATSHELL`, `MATMFFD`, `MatCreateShell()`, `MATTRANSPOSEVIRTUAL`, `MATSCHURCOMPLEMENT` 2005b77ba244SStefano Zampini @*/ 20069371c9d4SSatish Balay PetscErrorCode MatIsShell(Mat mat, PetscBool *flg) { 2007b77ba244SStefano Zampini PetscFunctionBegin; 2008b77ba244SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2009dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg, 2); 2010b77ba244SStefano Zampini *flg = (PetscBool)(mat->ops->destroy == MatDestroy_Shell); 2011b77ba244SStefano Zampini PetscFunctionReturn(0); 2012b77ba244SStefano Zampini } 2013