1be1d678aSKris Buschelman 2e51e0e81SBarry Smith /* 320563c6bSBarry Smith This provides a simple shell for Fortran (and C programmers) to 420563c6bSBarry Smith create a very simple matrix class for use with KSP without coding 5ed3cc1f0SBarry Smith much of anything. 6e51e0e81SBarry Smith */ 7e51e0e81SBarry Smith 8af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 945960306SStefano Zampini #include <petsc/private/vecimpl.h> 10e51e0e81SBarry Smith 1128f357e3SAlex Fikl struct _MatShellOps { 12976e8c5aSLisandro Dalcin /* 3 */ PetscErrorCode (*mult)(Mat,Vec,Vec); 13976e8c5aSLisandro Dalcin /* 5 */ PetscErrorCode (*multtranspose)(Mat,Vec,Vec); 14976e8c5aSLisandro Dalcin /* 17 */ PetscErrorCode (*getdiagonal)(Mat,Vec); 15976e8c5aSLisandro Dalcin /* 43 */ PetscErrorCode (*copy)(Mat,Mat,MatStructure); 16976e8c5aSLisandro Dalcin /* 60 */ PetscErrorCode (*destroy)(Mat); 1728f357e3SAlex Fikl }; 1828f357e3SAlex Fikl 19b77ba244SStefano Zampini struct _n_MatShellMatFunctionList { 20b77ba244SStefano Zampini PetscErrorCode (*symbolic)(Mat,Mat,Mat,void**); 21b77ba244SStefano Zampini PetscErrorCode (*numeric)(Mat,Mat,Mat,void*); 22b77ba244SStefano Zampini PetscErrorCode (*destroy)(void*); 23b77ba244SStefano Zampini MatProductType ptype; 24b77ba244SStefano Zampini char *composedname; /* string to identify routine with double dispatch */ 25b77ba244SStefano Zampini char *resultname; /* result matrix type */ 26b77ba244SStefano Zampini 27b77ba244SStefano Zampini struct _n_MatShellMatFunctionList *next; 28b77ba244SStefano Zampini }; 29b77ba244SStefano Zampini typedef struct _n_MatShellMatFunctionList *MatShellMatFunctionList; 30b77ba244SStefano Zampini 3128f357e3SAlex Fikl typedef struct { 3228f357e3SAlex Fikl struct _MatShellOps ops[1]; 332205254eSKarl Rupp 34b77ba244SStefano Zampini /* The user will manage the scaling and shifts for the MATSHELL, not the default */ 35b77ba244SStefano Zampini PetscBool managescalingshifts; 36b77ba244SStefano Zampini 37b77ba244SStefano Zampini /* support for MatScale, MatShift and MatMultAdd */ 38ef66eb69SBarry Smith PetscScalar vscale,vshift; 398fe8eb27SJed Brown Vec dshift; 408fe8eb27SJed Brown Vec left,right; 418fe8eb27SJed Brown Vec left_work,right_work; 425edf6516SJed Brown Vec left_add_work,right_add_work; 43b77ba244SStefano Zampini 44ef5c7bd2SStefano Zampini /* support for MatAXPY */ 459f137db4SBarry Smith Mat axpy; 469f137db4SBarry Smith PetscScalar axpy_vscale; 47b77ba244SStefano Zampini Vec axpy_left,axpy_right; 48ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 49b77ba244SStefano Zampini 5045960306SStefano Zampini /* support for ZeroRows/Columns operations */ 5145960306SStefano Zampini IS zrows; 5245960306SStefano Zampini IS zcols; 5345960306SStefano Zampini Vec zvals; 5445960306SStefano Zampini Vec zvals_w; 5545960306SStefano Zampini VecScatter zvals_sct_r; 5645960306SStefano Zampini VecScatter zvals_sct_c; 57b77ba244SStefano Zampini 58b77ba244SStefano Zampini /* MatMat operations */ 59b77ba244SStefano Zampini MatShellMatFunctionList matmat; 60b77ba244SStefano Zampini 61b77ba244SStefano Zampini /* user defined context */ 6220563c6bSBarry Smith void *ctx; 6388cf3e7dSBarry Smith } Mat_Shell; 640c0fd78eSBarry Smith 6545960306SStefano Zampini /* 6645960306SStefano Zampini Store and scale values on zeroed rows 6745960306SStefano Zampini xx = [x_1, 0], 0 on zeroed columns 6845960306SStefano Zampini */ 6945960306SStefano Zampini static PetscErrorCode MatShellPreZeroRight(Mat A,Vec x,Vec *xx) 7045960306SStefano Zampini { 7145960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 7245960306SStefano Zampini 7345960306SStefano Zampini PetscFunctionBegin; 7445960306SStefano Zampini *xx = x; 7545960306SStefano Zampini if (shell->zrows) { 765f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(shell->zvals_w,0.0)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals)); 8045960306SStefano Zampini } 8145960306SStefano Zampini if (shell->zcols) { 8245960306SStefano Zampini if (!shell->right_work) { 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&shell->right_work,NULL)); 8445960306SStefano Zampini } 855f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,shell->right_work)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecISSet(shell->right_work,shell->zcols,0.0)); 8745960306SStefano Zampini *xx = shell->right_work; 8845960306SStefano Zampini } 8945960306SStefano Zampini PetscFunctionReturn(0); 9045960306SStefano Zampini } 9145960306SStefano Zampini 9245960306SStefano Zampini /* Insert properly diagonally scaled values stored in MatShellPreZeroRight */ 9345960306SStefano Zampini static PetscErrorCode MatShellPostZeroLeft(Mat A,Vec x) 9445960306SStefano Zampini { 9545960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 9645960306SStefano Zampini 9745960306SStefano Zampini PetscFunctionBegin; 9845960306SStefano Zampini if (shell->zrows) { 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE)); 10145960306SStefano Zampini } 10245960306SStefano Zampini PetscFunctionReturn(0); 10345960306SStefano Zampini } 10445960306SStefano Zampini 10545960306SStefano Zampini /* 10645960306SStefano Zampini Store and scale values on zeroed rows 10745960306SStefano Zampini xx = [x_1, 0], 0 on zeroed rows 10845960306SStefano Zampini */ 10945960306SStefano Zampini static PetscErrorCode MatShellPreZeroLeft(Mat A,Vec x,Vec *xx) 11045960306SStefano Zampini { 11145960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 11245960306SStefano Zampini 11345960306SStefano Zampini PetscFunctionBegin; 11445960306SStefano Zampini *xx = NULL; 11545960306SStefano Zampini if (!shell->zrows) { 11645960306SStefano Zampini *xx = x; 11745960306SStefano Zampini } else { 11845960306SStefano Zampini if (!shell->left_work) { 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,NULL,&shell->left_work)); 12045960306SStefano Zampini } 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,shell->left_work)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(shell->zvals_w,0.0)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals)); 12845960306SStefano Zampini *xx = shell->left_work; 12945960306SStefano Zampini } 13045960306SStefano Zampini PetscFunctionReturn(0); 13145960306SStefano Zampini } 13245960306SStefano Zampini 13345960306SStefano Zampini /* Zero zero-columns contributions, sum contributions from properly scaled values stored in MatShellPreZeroLeft */ 13445960306SStefano Zampini static PetscErrorCode MatShellPostZeroRight(Mat A,Vec x) 13545960306SStefano Zampini { 13645960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 13745960306SStefano Zampini 13845960306SStefano Zampini PetscFunctionBegin; 13945960306SStefano Zampini if (shell->zcols) { 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecISSet(x,shell->zcols,0.0)); 14145960306SStefano Zampini } 14245960306SStefano Zampini if (shell->zrows) { 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE)); 14545960306SStefano Zampini } 14645960306SStefano Zampini PetscFunctionReturn(0); 14745960306SStefano Zampini } 14845960306SStefano Zampini 1498fe8eb27SJed Brown /* 1500c0fd78eSBarry Smith xx = diag(left)*x 1518fe8eb27SJed Brown */ 1528fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 1538fe8eb27SJed Brown { 1548fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1558fe8eb27SJed Brown 1568fe8eb27SJed Brown PetscFunctionBegin; 1570298fd71SBarry Smith *xx = NULL; 1588fe8eb27SJed Brown if (!shell->left) { 1598fe8eb27SJed Brown *xx = x; 1608fe8eb27SJed Brown } else { 1615f80ce2aSJacob Faibussowitsch if (!shell->left_work) CHKERRQ(VecDuplicate(shell->left,&shell->left_work)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->left_work,x,shell->left)); 1638fe8eb27SJed Brown *xx = shell->left_work; 1648fe8eb27SJed Brown } 1658fe8eb27SJed Brown PetscFunctionReturn(0); 1668fe8eb27SJed Brown } 1678fe8eb27SJed Brown 1680c0fd78eSBarry Smith /* 1690c0fd78eSBarry Smith xx = diag(right)*x 1700c0fd78eSBarry Smith */ 1718fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 1728fe8eb27SJed Brown { 1738fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1748fe8eb27SJed Brown 1758fe8eb27SJed Brown PetscFunctionBegin; 1760298fd71SBarry Smith *xx = NULL; 1778fe8eb27SJed Brown if (!shell->right) { 1788fe8eb27SJed Brown *xx = x; 1798fe8eb27SJed Brown } else { 1805f80ce2aSJacob Faibussowitsch if (!shell->right_work) CHKERRQ(VecDuplicate(shell->right,&shell->right_work)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->right_work,x,shell->right)); 1828fe8eb27SJed Brown *xx = shell->right_work; 1838fe8eb27SJed Brown } 1848fe8eb27SJed Brown PetscFunctionReturn(0); 1858fe8eb27SJed Brown } 1868fe8eb27SJed Brown 1870c0fd78eSBarry Smith /* 1880c0fd78eSBarry Smith x = diag(left)*x 1890c0fd78eSBarry Smith */ 1908fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 1918fe8eb27SJed Brown { 1928fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1938fe8eb27SJed Brown 1948fe8eb27SJed Brown PetscFunctionBegin; 1955f80ce2aSJacob Faibussowitsch if (shell->left) CHKERRQ(VecPointwiseMult(x,x,shell->left)); 1968fe8eb27SJed Brown PetscFunctionReturn(0); 1978fe8eb27SJed Brown } 1988fe8eb27SJed Brown 1990c0fd78eSBarry Smith /* 2000c0fd78eSBarry Smith x = diag(right)*x 2010c0fd78eSBarry Smith */ 2028fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 2038fe8eb27SJed Brown { 2048fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2058fe8eb27SJed Brown 2068fe8eb27SJed Brown PetscFunctionBegin; 2075f80ce2aSJacob Faibussowitsch if (shell->right) CHKERRQ(VecPointwiseMult(x,x,shell->right)); 2088fe8eb27SJed Brown PetscFunctionReturn(0); 2098fe8eb27SJed Brown } 2108fe8eb27SJed Brown 2110c0fd78eSBarry Smith /* 2120c0fd78eSBarry Smith Y = vscale*Y + diag(dshift)*X + vshift*X 2130c0fd78eSBarry Smith 2140c0fd78eSBarry Smith On input Y already contains A*x 2150c0fd78eSBarry Smith */ 2168fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 2178fe8eb27SJed Brown { 2188fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2198fe8eb27SJed Brown 2208fe8eb27SJed Brown PetscFunctionBegin; 2218fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 2228fe8eb27SJed Brown PetscInt i,m; 2238fe8eb27SJed Brown const PetscScalar *x,*d; 2248fe8eb27SJed Brown PetscScalar *y; 2255f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(X,&m)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(shell->dshift,&d)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(Y,&y)); 2298fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 2305f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(shell->dshift,&d)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(Y,&y)); 2330c0fd78eSBarry Smith } else { 2345f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(Y,shell->vscale)); 2358fe8eb27SJed Brown } 2365f80ce2aSJacob Faibussowitsch if (shell->vshift != 0.0) CHKERRQ(VecAXPY(Y,shell->vshift,X)); /* if test is for non-square matrices */ 2378fe8eb27SJed Brown PetscFunctionReturn(0); 2388fe8eb27SJed Brown } 239e51e0e81SBarry Smith 240789d8953SBarry Smith PetscErrorCode MatShellGetContext_Shell(Mat mat,void *ctx) 241789d8953SBarry Smith { 242789d8953SBarry Smith PetscFunctionBegin; 243789d8953SBarry Smith *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 244789d8953SBarry Smith PetscFunctionReturn(0); 245789d8953SBarry Smith } 246789d8953SBarry Smith 2479d225801SJed Brown /*@ 248a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 249b4fd4287SBarry Smith 250c7fcc2eaSBarry Smith Not Collective 251c7fcc2eaSBarry Smith 252b4fd4287SBarry Smith Input Parameter: 253b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 254b4fd4287SBarry Smith 255b4fd4287SBarry Smith Output Parameter: 256b4fd4287SBarry Smith . ctx - the user provided context 257b4fd4287SBarry Smith 25815091d37SBarry Smith Level: advanced 25915091d37SBarry Smith 26095452b02SPatrick Sanan Fortran Notes: 26195452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 262daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 263a62d957aSLois Curfman McInnes 264ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 2650bc0a719SBarry Smith @*/ 2668fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 267b4fd4287SBarry Smith { 2683a40ed3dSBarry Smith PetscFunctionBegin; 2690700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2704482741eSBarry Smith PetscValidPointer(ctx,2); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscUseMethod(mat,"MatShellGetContext_C",(Mat,void*),(mat,ctx))); 2723a40ed3dSBarry Smith PetscFunctionReturn(0); 273b4fd4287SBarry Smith } 274b4fd4287SBarry Smith 27545960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Local_Shell(Mat mat,PetscInt nr,PetscInt rows[],PetscInt nc,PetscInt cols[],PetscScalar diag,PetscBool rc) 27645960306SStefano Zampini { 27745960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)mat->data; 27845960306SStefano Zampini Vec x = NULL,b = NULL; 27945960306SStefano Zampini IS is1, is2; 28045960306SStefano Zampini const PetscInt *ridxs; 28145960306SStefano Zampini PetscInt *idxs,*gidxs; 28245960306SStefano Zampini PetscInt cum,rst,cst,i; 28345960306SStefano Zampini 28445960306SStefano Zampini PetscFunctionBegin; 28545960306SStefano Zampini if (!shell->zvals) { 2865f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(mat,NULL,&shell->zvals)); 28745960306SStefano Zampini } 28845960306SStefano Zampini if (!shell->zvals_w) { 2895f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(shell->zvals,&shell->zvals_w)); 29045960306SStefano Zampini } 2915f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(mat,&rst,NULL)); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRangeColumn(mat,&cst,NULL)); 29345960306SStefano Zampini 29445960306SStefano Zampini /* Expand/create index set of zeroed rows */ 2955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nr,&idxs)); 29645960306SStefano Zampini for (i = 0; i < nr; i++) idxs[i] = rows[i] + rst; 2975f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,nr,idxs,PETSC_OWN_POINTER,&is1)); 2985f80ce2aSJacob Faibussowitsch CHKERRQ(ISSort(is1)); 2995f80ce2aSJacob Faibussowitsch CHKERRQ(VecISSet(shell->zvals,is1,diag)); 30045960306SStefano Zampini if (shell->zrows) { 3015f80ce2aSJacob Faibussowitsch CHKERRQ(ISSum(shell->zrows,is1,&is2)); 3025f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&shell->zrows)); 3035f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is1)); 30445960306SStefano Zampini shell->zrows = is2; 30545960306SStefano Zampini } else shell->zrows = is1; 30645960306SStefano Zampini 30745960306SStefano Zampini /* Create scatters for diagonal values communications */ 3085f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&shell->zvals_sct_c)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&shell->zvals_sct_r)); 31045960306SStefano Zampini 31145960306SStefano Zampini /* row scatter: from/to left vector */ 3125f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(mat,&x,&b)); 3135f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(b,shell->zrows,shell->zvals_w,shell->zrows,&shell->zvals_sct_r)); 31445960306SStefano Zampini 31545960306SStefano Zampini /* col scatter: from right vector to left vector */ 3165f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(shell->zrows,&ridxs)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(shell->zrows,&nr)); 3185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nr,&gidxs)); 31945960306SStefano Zampini for (i = 0, cum = 0; i < nr; i++) { 32045960306SStefano Zampini if (ridxs[i] >= mat->cmap->N) continue; 32145960306SStefano Zampini gidxs[cum] = ridxs[i]; 32245960306SStefano Zampini cum++; 32345960306SStefano Zampini } 3245f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PetscObjectComm((PetscObject)mat),cum,gidxs,PETSC_OWN_POINTER,&is1)); 3255f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(x,is1,shell->zvals_w,is1,&shell->zvals_sct_c)); 3265f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is1)); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 32945960306SStefano Zampini 33045960306SStefano Zampini /* Expand/create index set of zeroed columns */ 33145960306SStefano Zampini if (rc) { 3325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nc,&idxs)); 33345960306SStefano Zampini for (i = 0; i < nc; i++) idxs[i] = cols[i] + cst; 3345f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,nc,idxs,PETSC_OWN_POINTER,&is1)); 3355f80ce2aSJacob Faibussowitsch CHKERRQ(ISSort(is1)); 33645960306SStefano Zampini if (shell->zcols) { 3375f80ce2aSJacob Faibussowitsch CHKERRQ(ISSum(shell->zcols,is1,&is2)); 3385f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&shell->zcols)); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is1)); 34045960306SStefano Zampini shell->zcols = is2; 34145960306SStefano Zampini } else shell->zcols = is1; 34245960306SStefano Zampini } 34345960306SStefano Zampini PetscFunctionReturn(0); 34445960306SStefano Zampini } 34545960306SStefano Zampini 34645960306SStefano Zampini static PetscErrorCode MatZeroRows_Shell(Mat mat,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 34745960306SStefano Zampini { 348ef5c7bd2SStefano Zampini Mat_Shell *shell = (Mat_Shell*)mat->data; 34945960306SStefano Zampini PetscInt nr, *lrows; 35045960306SStefano Zampini 35145960306SStefano Zampini PetscFunctionBegin; 35245960306SStefano Zampini if (x && b) { 35345960306SStefano Zampini Vec xt; 35445960306SStefano Zampini PetscScalar *vals; 35545960306SStefano Zampini PetscInt *gcols,i,st,nl,nc; 35645960306SStefano Zampini 3575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&gcols)); 35845960306SStefano Zampini for (i = 0, nc = 0; i < n; i++) if (rows[i] < mat->cmap->N) gcols[nc++] = rows[i]; 35945960306SStefano Zampini 3605f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(mat,&xt,NULL)); 3615f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,xt)); 3625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(nc,&vals)); 3635f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(xt,nc,gcols,vals,INSERT_VALUES)); /* xt = [x1, 0] */ 3645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(vals)); 3655f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(xt)); 3665f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(xt)); 3675f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(xt,-1.0,x)); /* xt = [0, x2] */ 36845960306SStefano Zampini 3695f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(xt,&st,NULL)); 3705f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(xt,&nl)); 3715f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xt,&vals)); 37245960306SStefano Zampini for (i = 0; i < nl; i++) { 37345960306SStefano Zampini PetscInt g = i + st; 37445960306SStefano Zampini if (g > mat->rmap->N) continue; 37545960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 3765f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValue(b,g,diag*vals[i],INSERT_VALUES)); 37745960306SStefano Zampini } 3785f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xt,&vals)); 3795f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(b)); 3805f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(b)); /* b = [b1, x2 * diag] */ 3815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xt)); 3825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(gcols)); 38345960306SStefano Zampini } 3845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutMapLocal(mat->rmap,n,rows,&nr,&lrows,NULL)); 3855f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsColumns_Local_Shell(mat,nr,lrows,0,NULL,diag,PETSC_FALSE)); 386ef5c7bd2SStefano Zampini if (shell->axpy) { 3875f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRows(shell->axpy,n,rows,0.0,NULL,NULL)); 388ef5c7bd2SStefano Zampini } 3895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(lrows)); 39045960306SStefano Zampini PetscFunctionReturn(0); 39145960306SStefano Zampini } 39245960306SStefano Zampini 39345960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Shell(Mat mat,PetscInt n,const PetscInt rowscols[],PetscScalar diag,Vec x,Vec b) 39445960306SStefano Zampini { 395ef5c7bd2SStefano Zampini Mat_Shell *shell = (Mat_Shell*)mat->data; 39645960306SStefano Zampini PetscInt *lrows, *lcols; 39745960306SStefano Zampini PetscInt nr, nc; 39845960306SStefano Zampini PetscBool congruent; 39945960306SStefano Zampini 40045960306SStefano Zampini PetscFunctionBegin; 40145960306SStefano Zampini if (x && b) { 40245960306SStefano Zampini Vec xt, bt; 40345960306SStefano Zampini PetscScalar *vals; 40445960306SStefano Zampini PetscInt *grows,*gcols,i,st,nl; 40545960306SStefano Zampini 4065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(n,&grows,n,&gcols)); 40745960306SStefano Zampini for (i = 0, nr = 0; i < n; i++) if (rowscols[i] < mat->rmap->N) grows[nr++] = rowscols[i]; 40845960306SStefano Zampini for (i = 0, nc = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) gcols[nc++] = rowscols[i]; 4095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(n,&vals)); 41045960306SStefano Zampini 4115f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(mat,&xt,&bt)); 4125f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,xt)); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(xt,nc,gcols,vals,INSERT_VALUES)); /* xt = [x1, 0] */ 4145f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(xt)); 4155f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(xt)); 4165f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(xt,-1.0,x)); /* xt = [0, -x2] */ 4175f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(mat,xt,bt)); /* bt = [-A12*x2,-A22*x2] */ 4185f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(bt,nr,grows,vals,INSERT_VALUES)); /* bt = [-A12*x2,0] */ 4195f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(bt)); 4205f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(bt)); 4215f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b,1.0,bt)); /* b = [b1 - A12*x2, b2] */ 4225f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(bt,nr,grows,vals,INSERT_VALUES)); /* b = [b1 - A12*x2, 0] */ 4235f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(bt)); 4245f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(bt)); 4255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(vals)); 42645960306SStefano Zampini 4275f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(xt,&st,NULL)); 4285f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(xt,&nl)); 4295f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xt,&vals)); 43045960306SStefano Zampini for (i = 0; i < nl; i++) { 43145960306SStefano Zampini PetscInt g = i + st; 43245960306SStefano Zampini if (g > mat->rmap->N) continue; 43345960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 4345f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValue(b,g,-diag*vals[i],INSERT_VALUES)); 43545960306SStefano Zampini } 4365f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xt,&vals)); 4375f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(b)); 4385f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(b)); /* b = [b1 - A12*x2, x2 * diag] */ 4395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xt)); 4405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&bt)); 4415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(grows,gcols)); 44245960306SStefano Zampini } 4435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutMapLocal(mat->rmap,n,rowscols,&nr,&lrows,NULL)); 4445f80ce2aSJacob Faibussowitsch CHKERRQ(MatHasCongruentLayouts(mat,&congruent)); 44545960306SStefano Zampini if (congruent) { 44645960306SStefano Zampini nc = nr; 44745960306SStefano Zampini lcols = lrows; 44845960306SStefano Zampini } else { /* MatZeroRowsColumns implicitly assumes the rowscols indices are for a square matrix, here we handle a more general case */ 44945960306SStefano Zampini PetscInt i,nt,*t; 45045960306SStefano Zampini 4515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&t)); 45245960306SStefano Zampini for (i = 0, nt = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) t[nt++] = rowscols[i]; 4535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutMapLocal(mat->cmap,nt,t,&nc,&lcols,NULL)); 4545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(t)); 45545960306SStefano Zampini } 4565f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsColumns_Local_Shell(mat,nr,lrows,nc,lcols,diag,PETSC_TRUE)); 45745960306SStefano Zampini if (!congruent) { 4585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(lcols)); 45945960306SStefano Zampini } 4605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(lrows)); 461ef5c7bd2SStefano Zampini if (shell->axpy) { 4625f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsColumns(shell->axpy,n,rowscols,0.0,NULL,NULL)); 463ef5c7bd2SStefano Zampini } 46445960306SStefano Zampini PetscFunctionReturn(0); 46545960306SStefano Zampini } 46645960306SStefano Zampini 467dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 468e51e0e81SBarry Smith { 469bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 470b77ba244SStefano Zampini MatShellMatFunctionList matmat; 471ed3cc1f0SBarry Smith 4723a40ed3dSBarry Smith PetscFunctionBegin; 47328f357e3SAlex Fikl if (shell->ops->destroy) { 4745f80ce2aSJacob Faibussowitsch CHKERRQ((*shell->ops->destroy)(mat)); 475bf0cc555SLisandro Dalcin } 4765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(shell->ops,sizeof(struct _MatShellOps))); 4775f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->left)); 4785f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->right)); 4795f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->dshift)); 4805f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->left_work)); 4815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->right_work)); 4825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->left_add_work)); 4835f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->right_add_work)); 4845f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->axpy_left)); 4855f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->axpy_right)); 4865f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&shell->axpy)); 4875f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->zvals_w)); 4885f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->zvals)); 4895f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&shell->zvals_sct_c)); 4905f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&shell->zvals_sct_r)); 4915f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&shell->zrows)); 4925f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&shell->zcols)); 493b77ba244SStefano Zampini 494b77ba244SStefano Zampini matmat = shell->matmat; 495b77ba244SStefano Zampini while (matmat) { 496b77ba244SStefano Zampini MatShellMatFunctionList next = matmat->next; 497b77ba244SStefano Zampini 4985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,matmat->composedname,NULL)); 4995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(matmat->composedname)); 5005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(matmat->resultname)); 5015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(matmat)); 502b77ba244SStefano Zampini matmat = next; 503b77ba244SStefano Zampini } 5045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatShellGetContext_C",NULL)); 5055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetContext_C",NULL)); 5065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetVecType_C",NULL)); 5075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetManageScalingShifts_C",NULL)); 5085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetOperation_C",NULL)); 5095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatShellGetOperation_C",NULL)); 5105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetMatProductOperation_C",NULL)); 5115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(mat->data)); 512b77ba244SStefano Zampini PetscFunctionReturn(0); 513b77ba244SStefano Zampini } 514b77ba244SStefano Zampini 515b77ba244SStefano Zampini typedef struct { 516b77ba244SStefano Zampini PetscErrorCode (*numeric)(Mat,Mat,Mat,void*); 517b77ba244SStefano Zampini PetscErrorCode (*destroy)(void*); 518b77ba244SStefano Zampini void *userdata; 519b77ba244SStefano Zampini Mat B; 520b77ba244SStefano Zampini Mat Bt; 521b77ba244SStefano Zampini Mat axpy; 522b77ba244SStefano Zampini } MatMatDataShell; 523b77ba244SStefano Zampini 524b77ba244SStefano Zampini static PetscErrorCode DestroyMatMatDataShell(void *data) 525b77ba244SStefano Zampini { 526b77ba244SStefano Zampini MatMatDataShell *mmdata = (MatMatDataShell *)data; 527b77ba244SStefano Zampini 528b77ba244SStefano Zampini PetscFunctionBegin; 529b77ba244SStefano Zampini if (mmdata->destroy) { 5305f80ce2aSJacob Faibussowitsch CHKERRQ((*mmdata->destroy)(mmdata->userdata)); 531b77ba244SStefano Zampini } 5325f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mmdata->B)); 5335f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mmdata->Bt)); 5345f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mmdata->axpy)); 5355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(mmdata)); 536b77ba244SStefano Zampini PetscFunctionReturn(0); 537b77ba244SStefano Zampini } 538b77ba244SStefano Zampini 539b77ba244SStefano Zampini static PetscErrorCode MatProductNumeric_Shell_X(Mat D) 540b77ba244SStefano Zampini { 541b77ba244SStefano Zampini Mat_Product *product; 542b77ba244SStefano Zampini Mat A, B; 543b77ba244SStefano Zampini MatMatDataShell *mdata; 544b77ba244SStefano Zampini PetscScalar zero = 0.0; 545b77ba244SStefano Zampini 546b77ba244SStefano Zampini PetscFunctionBegin; 547b77ba244SStefano Zampini MatCheckProduct(D,1); 548b77ba244SStefano Zampini product = D->product; 549*28b400f6SJacob Faibussowitsch PetscCheck(product->data,PetscObjectComm((PetscObject)D),PETSC_ERR_PLIB,"Product data empty"); 550b77ba244SStefano Zampini A = product->A; 551b77ba244SStefano Zampini B = product->B; 552b77ba244SStefano Zampini mdata = (MatMatDataShell*)product->data; 553b77ba244SStefano Zampini if (mdata->numeric) { 554b77ba244SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 555b77ba244SStefano Zampini PetscErrorCode (*stashsym)(Mat) = D->ops->productsymbolic; 556b77ba244SStefano Zampini PetscErrorCode (*stashnum)(Mat) = D->ops->productnumeric; 557b77ba244SStefano Zampini PetscBool useBmdata = PETSC_FALSE, newB = PETSC_TRUE; 558b77ba244SStefano Zampini 559b77ba244SStefano Zampini if (shell->managescalingshifts) { 5602c71b3e2SJacob Faibussowitsch PetscCheckFalse(shell->zcols || shell->zrows,PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProduct not supported with zeroed rows/columns"); 561b77ba244SStefano Zampini if (shell->right || shell->left) { 562b77ba244SStefano Zampini useBmdata = PETSC_TRUE; 563b77ba244SStefano Zampini if (!mdata->B) { 5645f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(B,MAT_SHARE_NONZERO_PATTERN,&mdata->B)); 565b77ba244SStefano Zampini } else { 566b77ba244SStefano Zampini newB = PETSC_FALSE; 567b77ba244SStefano Zampini } 5685f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(B,mdata->B,SAME_NONZERO_PATTERN)); 569b77ba244SStefano Zampini } 570b77ba244SStefano Zampini switch (product->type) { 571b77ba244SStefano Zampini case MATPRODUCT_AB: /* s L A R B + v L R B + L D R B */ 572b77ba244SStefano Zampini if (shell->right) { 5735f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(mdata->B,shell->right,NULL)); 574b77ba244SStefano Zampini } 575b77ba244SStefano Zampini break; 576b77ba244SStefano Zampini case MATPRODUCT_AtB: /* s R A^t L B + v R L B + R D L B */ 577b77ba244SStefano Zampini if (shell->left) { 5785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(mdata->B,shell->left,NULL)); 579b77ba244SStefano Zampini } 580b77ba244SStefano Zampini break; 581b77ba244SStefano Zampini case MATPRODUCT_ABt: /* s L A R B^t + v L R B^t + L D R B^t */ 582b77ba244SStefano Zampini if (shell->right) { 5835f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(mdata->B,NULL,shell->right)); 584b77ba244SStefano Zampini } 585b77ba244SStefano Zampini break; 586b77ba244SStefano Zampini case MATPRODUCT_RARt: /* s B L A R B^t + v B L R B^t + B L D R B^t */ 587b77ba244SStefano Zampini if (shell->right && shell->left) { 588b77ba244SStefano Zampini PetscBool flg; 589b77ba244SStefano Zampini 5905f80ce2aSJacob Faibussowitsch CHKERRQ(VecEqual(shell->right,shell->left,&flg)); 591*28b400f6SJacob Faibussowitsch 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,((PetscObject)B)->type_name); 592b77ba244SStefano Zampini } 593b77ba244SStefano Zampini if (shell->right) { 5945f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(mdata->B,NULL,shell->right)); 595b77ba244SStefano Zampini } 596b77ba244SStefano Zampini break; 597b77ba244SStefano Zampini case MATPRODUCT_PtAP: /* s B^t L A R B + v B^t L R B + B^t L D R B */ 598b77ba244SStefano Zampini if (shell->right && shell->left) { 599b77ba244SStefano Zampini PetscBool flg; 600b77ba244SStefano Zampini 6015f80ce2aSJacob Faibussowitsch CHKERRQ(VecEqual(shell->right,shell->left,&flg)); 602*28b400f6SJacob Faibussowitsch 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,((PetscObject)B)->type_name); 603b77ba244SStefano Zampini } 604b77ba244SStefano Zampini if (shell->right) { 6055f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(mdata->B,shell->right,NULL)); 606b77ba244SStefano Zampini } 607b77ba244SStefano Zampini break; 60898921bdaSJacob 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); 609b77ba244SStefano Zampini } 610b77ba244SStefano Zampini } 611b77ba244SStefano Zampini /* allow the user to call MatMat operations on D */ 612b77ba244SStefano Zampini D->product = NULL; 613b77ba244SStefano Zampini D->ops->productsymbolic = NULL; 614b77ba244SStefano Zampini D->ops->productnumeric = NULL; 615b77ba244SStefano Zampini 6165f80ce2aSJacob Faibussowitsch CHKERRQ((*mdata->numeric)(A,useBmdata ? mdata->B : B,D,mdata->userdata)); 617b77ba244SStefano Zampini 618b77ba244SStefano Zampini /* clear any leftover user data and restore D pointers */ 6195f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(D)); 620b77ba244SStefano Zampini D->ops->productsymbolic = stashsym; 621b77ba244SStefano Zampini D->ops->productnumeric = stashnum; 622b77ba244SStefano Zampini D->product = product; 623b77ba244SStefano Zampini 624b77ba244SStefano Zampini if (shell->managescalingshifts) { 6255f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(D,shell->vscale)); 626b77ba244SStefano Zampini switch (product->type) { 627b77ba244SStefano Zampini case MATPRODUCT_AB: /* s L A R B + v L R B + L D R B */ 628b77ba244SStefano Zampini case MATPRODUCT_ABt: /* s L A R B^t + v L R B^t + L D R B^t */ 629b77ba244SStefano Zampini if (shell->left) { 6305f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(D,shell->left,NULL)); 631b77ba244SStefano Zampini if (shell->dshift || shell->vshift != zero) { 6325f80ce2aSJacob Faibussowitsch if (!shell->left_work) CHKERRQ(MatCreateVecs(A,NULL,&shell->left_work)); 633b77ba244SStefano Zampini if (shell->dshift) { 6345f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(shell->dshift,shell->left_work)); 6355f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(shell->left_work,shell->vshift)); 6365f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->left_work,shell->left_work,shell->left)); 637b77ba244SStefano Zampini } else { 6385f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(shell->left_work,shell->vshift)); 639b77ba244SStefano Zampini } 640b77ba244SStefano Zampini if (product->type == MATPRODUCT_ABt) { 641b77ba244SStefano Zampini MatReuse reuse = mdata->Bt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX; 642b77ba244SStefano Zampini MatStructure str = mdata->Bt ? SUBSET_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 643b77ba244SStefano Zampini 6445f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(mdata->B,reuse,&mdata->Bt)); 6455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(mdata->Bt,shell->left_work,NULL)); 6465f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(D,1.0,mdata->Bt,str)); 647b77ba244SStefano Zampini } else { 648b77ba244SStefano Zampini MatStructure str = newB ? DIFFERENT_NONZERO_PATTERN : SUBSET_NONZERO_PATTERN; 649b77ba244SStefano Zampini 6505f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(mdata->B,shell->left_work,NULL)); 6515f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(D,1.0,mdata->B,str)); 652b77ba244SStefano Zampini } 653b77ba244SStefano Zampini } 654b77ba244SStefano Zampini } 655b77ba244SStefano Zampini break; 656b77ba244SStefano Zampini case MATPRODUCT_AtB: /* s R A^t L B + v R L B + R D L B */ 657b77ba244SStefano Zampini if (shell->right) { 6585f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(D,shell->right,NULL)); 659b77ba244SStefano Zampini if (shell->dshift || shell->vshift != zero) { 660b77ba244SStefano Zampini MatStructure str = newB ? DIFFERENT_NONZERO_PATTERN : SUBSET_NONZERO_PATTERN; 661b77ba244SStefano Zampini 6625f80ce2aSJacob Faibussowitsch if (!shell->right_work) CHKERRQ(MatCreateVecs(A,&shell->right_work,NULL)); 663b77ba244SStefano Zampini if (shell->dshift) { 6645f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(shell->dshift,shell->right_work)); 6655f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(shell->right_work,shell->vshift)); 6665f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->right_work,shell->right_work,shell->right)); 667b77ba244SStefano Zampini } else { 6685f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(shell->right_work,shell->vshift)); 669b77ba244SStefano Zampini } 6705f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(mdata->B,shell->right_work,NULL)); 6715f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(D,1.0,mdata->B,str)); 672b77ba244SStefano Zampini } 673b77ba244SStefano Zampini } 674b77ba244SStefano Zampini break; 675b77ba244SStefano Zampini case MATPRODUCT_PtAP: /* s B^t L A R B + v B^t L R B + B^t L D R B */ 676b77ba244SStefano Zampini case MATPRODUCT_RARt: /* s B L A R B^t + v B L R B^t + B L D R B^t */ 6772c71b3e2SJacob Faibussowitsch PetscCheckFalse(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,((PetscObject)B)->type_name); 678b77ba244SStefano Zampini break; 67998921bdaSJacob 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); 680b77ba244SStefano Zampini } 681b77ba244SStefano Zampini if (shell->axpy && shell->axpy_vscale != zero) { 682b77ba244SStefano Zampini Mat X; 683b77ba244SStefano Zampini PetscObjectState axpy_state; 684b77ba244SStefano Zampini MatStructure str = DIFFERENT_NONZERO_PATTERN; /* not sure it is safe to ever use SUBSET_NONZERO_PATTERN */ 685b77ba244SStefano Zampini 6865f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(shell->axpy,&X)); 6875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)X,&axpy_state)); 6882c71b3e2SJacob Faibussowitsch PetscCheckFalse(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,...)"); 689b77ba244SStefano Zampini if (!mdata->axpy) { 690b77ba244SStefano Zampini str = DIFFERENT_NONZERO_PATTERN; 6915f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreate(shell->axpy,B,NULL,&mdata->axpy)); 6925f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(mdata->axpy,product->type)); 6935f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(mdata->axpy)); 6945f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(mdata->axpy)); 695b77ba244SStefano Zampini } else { /* May be that shell->axpy has changed */ 696b77ba244SStefano Zampini PetscBool flg; 697b77ba244SStefano Zampini 6985f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductReplaceMats(shell->axpy,B,NULL,mdata->axpy)); 6995f80ce2aSJacob Faibussowitsch CHKERRQ(MatHasOperation(mdata->axpy,MATOP_PRODUCTSYMBOLIC,&flg)); 700b77ba244SStefano Zampini if (!flg) { 701b77ba244SStefano Zampini str = DIFFERENT_NONZERO_PATTERN; 7025f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(mdata->axpy)); 7035f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(mdata->axpy)); 704b77ba244SStefano Zampini } 705b77ba244SStefano Zampini } 7065f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(mdata->axpy)); 7075f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(D,shell->axpy_vscale,mdata->axpy,str)); 708b77ba244SStefano Zampini } 709b77ba244SStefano Zampini } 710b77ba244SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)D),PETSC_ERR_PLIB,"Missing numeric operation"); 711b77ba244SStefano Zampini PetscFunctionReturn(0); 712b77ba244SStefano Zampini } 713b77ba244SStefano Zampini 714b77ba244SStefano Zampini static PetscErrorCode MatProductSymbolic_Shell_X(Mat D) 715b77ba244SStefano Zampini { 716b77ba244SStefano Zampini Mat_Product *product; 717b77ba244SStefano Zampini Mat A,B; 718b77ba244SStefano Zampini MatShellMatFunctionList matmat; 719b77ba244SStefano Zampini Mat_Shell *shell; 720b77ba244SStefano Zampini PetscBool flg; 721b77ba244SStefano Zampini char composedname[256]; 722b77ba244SStefano Zampini MatMatDataShell *mdata; 723b77ba244SStefano Zampini 724b77ba244SStefano Zampini PetscFunctionBegin; 725b77ba244SStefano Zampini MatCheckProduct(D,1); 726b77ba244SStefano Zampini product = D->product; 727*28b400f6SJacob Faibussowitsch PetscCheck(!product->data,PetscObjectComm((PetscObject)D),PETSC_ERR_PLIB,"Product data not empty"); 728b77ba244SStefano Zampini A = product->A; 729b77ba244SStefano Zampini B = product->B; 730b77ba244SStefano Zampini shell = (Mat_Shell*)A->data; 731b77ba244SStefano Zampini matmat = shell->matmat; 7325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(composedname,sizeof(composedname),"MatProductSetFromOptions_%s_%s_C",((PetscObject)A)->type_name,((PetscObject)B)->type_name)); 733b77ba244SStefano Zampini while (matmat) { 7345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(composedname,matmat->composedname,&flg)); 735b77ba244SStefano Zampini flg = (PetscBool)(flg && (matmat->ptype == product->type)); 736b77ba244SStefano Zampini if (flg) break; 737b77ba244SStefano Zampini matmat = matmat->next; 738b77ba244SStefano Zampini } 739*28b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Composedname \"%s\" for product type %s not found",composedname,MatProductTypes[product->type]); 740b77ba244SStefano Zampini switch (product->type) { 741b77ba244SStefano Zampini case MATPRODUCT_AB: 7425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(D,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N)); 743b77ba244SStefano Zampini break; 744b77ba244SStefano Zampini case MATPRODUCT_AtB: 7455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(D,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N)); 746b77ba244SStefano Zampini break; 747b77ba244SStefano Zampini case MATPRODUCT_ABt: 7485f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(D,A->rmap->n,B->rmap->n,A->rmap->N,B->rmap->N)); 749b77ba244SStefano Zampini break; 750b77ba244SStefano Zampini case MATPRODUCT_RARt: 7515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(D,B->rmap->n,B->rmap->n,B->rmap->N,B->rmap->N)); 752b77ba244SStefano Zampini break; 753b77ba244SStefano Zampini case MATPRODUCT_PtAP: 7545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(D,B->cmap->n,B->cmap->n,B->cmap->N,B->cmap->N)); 755b77ba244SStefano Zampini break; 75698921bdaSJacob 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); 757b77ba244SStefano Zampini } 758b77ba244SStefano Zampini /* respect users who passed in a matrix for which resultname is the base type */ 759b77ba244SStefano Zampini if (matmat->resultname) { 7605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)D,matmat->resultname,&flg)); 761b77ba244SStefano Zampini if (!flg) { 7625f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(D,matmat->resultname)); 763b77ba244SStefano Zampini } 764b77ba244SStefano Zampini } 765b77ba244SStefano Zampini /* If matrix type was not set or different, we need to reset this pointers */ 766b77ba244SStefano Zampini D->ops->productsymbolic = MatProductSymbolic_Shell_X; 767b77ba244SStefano Zampini D->ops->productnumeric = MatProductNumeric_Shell_X; 768b77ba244SStefano Zampini /* attach product data */ 7695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&mdata)); 770b77ba244SStefano Zampini mdata->numeric = matmat->numeric; 771b77ba244SStefano Zampini mdata->destroy = matmat->destroy; 772b77ba244SStefano Zampini if (matmat->symbolic) { 7735f80ce2aSJacob Faibussowitsch CHKERRQ((*matmat->symbolic)(A,B,D,&mdata->userdata)); 774b77ba244SStefano Zampini } else { /* call general setup if symbolic operation not provided */ 7755f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(D)); 776b77ba244SStefano Zampini } 777*28b400f6SJacob Faibussowitsch PetscCheck(D->product,PetscObjectComm((PetscObject)D),PETSC_ERR_COR,"Product disappeared after user symbolic phase"); 778*28b400f6SJacob Faibussowitsch PetscCheck(!D->product->data,PetscObjectComm((PetscObject)D),PETSC_ERR_COR,"Product data not empty after user symbolic phase"); 779b77ba244SStefano Zampini D->product->data = mdata; 780b77ba244SStefano Zampini D->product->destroy = DestroyMatMatDataShell; 781b77ba244SStefano Zampini /* Be sure to reset these pointers if the user did something unexpected */ 782b77ba244SStefano Zampini D->ops->productsymbolic = MatProductSymbolic_Shell_X; 783b77ba244SStefano Zampini D->ops->productnumeric = MatProductNumeric_Shell_X; 784b77ba244SStefano Zampini PetscFunctionReturn(0); 785b77ba244SStefano Zampini } 786b77ba244SStefano Zampini 787b77ba244SStefano Zampini static PetscErrorCode MatProductSetFromOptions_Shell_X(Mat D) 788b77ba244SStefano Zampini { 789b77ba244SStefano Zampini Mat_Product *product; 790b77ba244SStefano Zampini Mat A,B; 791b77ba244SStefano Zampini MatShellMatFunctionList matmat; 792b77ba244SStefano Zampini Mat_Shell *shell; 793b77ba244SStefano Zampini PetscBool flg; 794b77ba244SStefano Zampini char composedname[256]; 795b77ba244SStefano Zampini 796b77ba244SStefano Zampini PetscFunctionBegin; 797b77ba244SStefano Zampini MatCheckProduct(D,1); 798b77ba244SStefano Zampini product = D->product; 799b77ba244SStefano Zampini A = product->A; 800b77ba244SStefano Zampini B = product->B; 8015f80ce2aSJacob Faibussowitsch CHKERRQ(MatIsShell(A,&flg)); 802b77ba244SStefano Zampini if (!flg) PetscFunctionReturn(0); 803b77ba244SStefano Zampini shell = (Mat_Shell*)A->data; 804b77ba244SStefano Zampini matmat = shell->matmat; 8055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(composedname,sizeof(composedname),"MatProductSetFromOptions_%s_%s_C",((PetscObject)A)->type_name,((PetscObject)B)->type_name)); 806b77ba244SStefano Zampini while (matmat) { 8075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(composedname,matmat->composedname,&flg)); 808b77ba244SStefano Zampini flg = (PetscBool)(flg && (matmat->ptype == product->type)); 809b77ba244SStefano Zampini if (flg) break; 810b77ba244SStefano Zampini matmat = matmat->next; 811b77ba244SStefano Zampini } 812b77ba244SStefano Zampini if (flg) { D->ops->productsymbolic = MatProductSymbolic_Shell_X; } 8135f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscInfo(D," symbolic product %s not registered for product type %s\n",composedname,MatProductTypes[product->type])); 814b77ba244SStefano Zampini PetscFunctionReturn(0); 815b77ba244SStefano Zampini } 816b77ba244SStefano Zampini 817b77ba244SStefano Zampini static PetscErrorCode MatShellSetMatProductOperation_Private(Mat A,MatProductType ptype,PetscErrorCode (*symbolic)(Mat,Mat,Mat,void**),PetscErrorCode (*numeric)(Mat,Mat,Mat,void*),PetscErrorCode (*destroy)(void*),char *composedname,const char *resultname) 818b77ba244SStefano Zampini { 819b77ba244SStefano Zampini PetscBool flg; 820b77ba244SStefano Zampini Mat_Shell *shell; 821b77ba244SStefano Zampini MatShellMatFunctionList matmat; 822b77ba244SStefano Zampini 823b77ba244SStefano Zampini PetscFunctionBegin; 824*28b400f6SJacob Faibussowitsch PetscCheck(numeric,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing numeric routine"); 825*28b400f6SJacob Faibussowitsch PetscCheck(composedname,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing composed name"); 826b77ba244SStefano Zampini 827b77ba244SStefano Zampini /* add product callback */ 828b77ba244SStefano Zampini shell = (Mat_Shell*)A->data; 829b77ba244SStefano Zampini matmat = shell->matmat; 830b77ba244SStefano Zampini if (!matmat) { 8315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&shell->matmat)); 832b77ba244SStefano Zampini matmat = shell->matmat; 833b77ba244SStefano Zampini } else { 834b77ba244SStefano Zampini MatShellMatFunctionList entry = matmat; 835b77ba244SStefano Zampini while (entry) { 8365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(composedname,entry->composedname,&flg)); 837b77ba244SStefano Zampini flg = (PetscBool)(flg && (entry->ptype == ptype)); 838843d480fSStefano Zampini if (flg) goto set; 839b77ba244SStefano Zampini matmat = entry; 840b77ba244SStefano Zampini entry = entry->next; 841b77ba244SStefano Zampini } 8425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&matmat->next)); 843b77ba244SStefano Zampini matmat = matmat->next; 844b77ba244SStefano Zampini } 845b77ba244SStefano Zampini 846843d480fSStefano Zampini set: 847b77ba244SStefano Zampini matmat->symbolic = symbolic; 848b77ba244SStefano Zampini matmat->numeric = numeric; 849b77ba244SStefano Zampini matmat->destroy = destroy; 850b77ba244SStefano Zampini matmat->ptype = ptype; 8515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(matmat->composedname)); 8525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(matmat->resultname)); 8535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(composedname,&matmat->composedname)); 8545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(resultname,&matmat->resultname)); 8555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"Composing %s for product type %s with result %s\n",matmat->composedname,MatProductTypes[matmat->ptype],matmat->resultname ? matmat->resultname : "not specified")); 8565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,matmat->composedname,MatProductSetFromOptions_Shell_X)); 857b77ba244SStefano Zampini PetscFunctionReturn(0); 858b77ba244SStefano Zampini } 859b77ba244SStefano Zampini 860b77ba244SStefano Zampini /*@C 861b77ba244SStefano Zampini MatShellSetMatProductOperation - Allows user to set a matrix matrix operation for a shell matrix. 862b77ba244SStefano Zampini 863b77ba244SStefano Zampini Logically Collective on Mat 864b77ba244SStefano Zampini 865b77ba244SStefano Zampini Input Parameters: 866b77ba244SStefano Zampini + A - the shell matrix 867b77ba244SStefano Zampini . ptype - the product type 868b77ba244SStefano Zampini . symbolic - the function for the symbolic phase (can be NULL) 869b77ba244SStefano Zampini . numeric - the function for the numerical phase 870b77ba244SStefano Zampini . destroy - the function for the destruction of the needed data generated during the symbolic phase (can be NULL) 871b77ba244SStefano Zampini . Btype - the matrix type for the matrix to be multiplied against 872b77ba244SStefano Zampini - Ctype - the matrix type for the result (can be NULL) 873b77ba244SStefano Zampini 874b77ba244SStefano Zampini Level: advanced 875b77ba244SStefano Zampini 876b77ba244SStefano Zampini Usage: 877b77ba244SStefano Zampini $ extern PetscErrorCode usersymbolic(Mat,Mat,Mat,void**); 878b77ba244SStefano Zampini $ extern PetscErrorCode usernumeric(Mat,Mat,Mat,void*); 879b77ba244SStefano Zampini $ extern PetscErrorCode userdestroy(void*); 880b77ba244SStefano Zampini $ MatCreateShell(comm,m,n,M,N,ctx,&A); 881b77ba244SStefano Zampini $ MatShellSetMatProductOperation(A,MATPRODUCT_AB,usersymbolic,usernumeric,userdestroy,MATSEQAIJ,MATDENSE); 882b77ba244SStefano Zampini $ [ create B of type SEQAIJ etc..] 883b77ba244SStefano Zampini $ MatProductCreate(A,B,NULL,&C); 884b77ba244SStefano Zampini $ MatProductSetType(C,MATPRODUCT_AB); 885b77ba244SStefano Zampini $ MatProductSetFromOptions(C); 886b77ba244SStefano Zampini $ MatProductSymbolic(C); -> actually runs the user defined symbolic operation 887b77ba244SStefano Zampini $ MatProductNumeric(C); -> actually runs the user defined numeric operation 888b77ba244SStefano Zampini $ [ use C = A*B ] 889b77ba244SStefano Zampini 890b77ba244SStefano Zampini Notes: 891b77ba244SStefano Zampini MATPRODUCT_ABC is not supported yet. Not supported in Fortran. 892b77ba244SStefano Zampini If the symbolic phase is not specified, MatSetUp() is called on the result matrix that must have its type set if Ctype is NULL. 893b77ba244SStefano Zampini Any additional data needed by the matrix product needs to be returned during the symbolic phase and destroyed with the destroy callback. 894b77ba244SStefano Zampini PETSc will take care of calling the user-defined callbacks. 895b77ba244SStefano Zampini It is allowed to specify the same callbacks for different Btype matrix types. 896b77ba244SStefano Zampini The couple (Btype,ptype) uniquely identifies the operation: the last specified callbacks takes precedence. 897b77ba244SStefano Zampini 898b77ba244SStefano Zampini .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatProductType, MatType, MatSetUp() 899b77ba244SStefano Zampini @*/ 900b77ba244SStefano Zampini PetscErrorCode MatShellSetMatProductOperation(Mat A,MatProductType ptype,PetscErrorCode (*symbolic)(Mat,Mat,Mat,void**),PetscErrorCode (*numeric)(Mat,Mat,Mat,void*),PetscErrorCode (*destroy)(void *),MatType Btype,MatType Ctype) 901b77ba244SStefano Zampini { 902b77ba244SStefano Zampini PetscFunctionBegin; 903b77ba244SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 904b77ba244SStefano Zampini PetscValidLogicalCollectiveEnum(A,ptype,2); 9052c71b3e2SJacob Faibussowitsch PetscCheckFalse(ptype == MATPRODUCT_ABC,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not for product type %s",MatProductTypes[ptype]); 906*28b400f6SJacob Faibussowitsch PetscCheck(numeric,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Missing numeric routine, argument 4"); 907b77ba244SStefano Zampini PetscValidPointer(Btype,6); 908b77ba244SStefano Zampini if (Ctype) PetscValidPointer(Ctype,7); 9095f80ce2aSJacob Faibussowitsch CHKERRQ(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))); 910b77ba244SStefano Zampini PetscFunctionReturn(0); 911b77ba244SStefano Zampini } 912b77ba244SStefano Zampini 913b77ba244SStefano Zampini PetscErrorCode MatShellSetMatProductOperation_Shell(Mat A,MatProductType ptype,PetscErrorCode (*symbolic)(Mat,Mat,Mat,void**),PetscErrorCode (*numeric)(Mat,Mat,Mat,void*),PetscErrorCode (*destroy)(void *),MatType Btype,MatType Ctype) 914b77ba244SStefano Zampini { 915b77ba244SStefano Zampini PetscBool flg; 916b77ba244SStefano Zampini char composedname[256]; 917b77ba244SStefano Zampini MatRootName Bnames = MatRootNameList, Cnames = MatRootNameList; 918b77ba244SStefano Zampini PetscMPIInt size; 919b77ba244SStefano Zampini 920b77ba244SStefano Zampini PetscFunctionBegin; 921b77ba244SStefano Zampini PetscValidType(A,1); 922b77ba244SStefano Zampini while (Bnames) { /* user passed in the root name */ 9235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(Btype,Bnames->rname,&flg)); 924b77ba244SStefano Zampini if (flg) break; 925b77ba244SStefano Zampini Bnames = Bnames->next; 926b77ba244SStefano Zampini } 927b77ba244SStefano Zampini while (Cnames) { /* user passed in the root name */ 9285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(Ctype,Cnames->rname,&flg)); 929b77ba244SStefano Zampini if (flg) break; 930b77ba244SStefano Zampini Cnames = Cnames->next; 931b77ba244SStefano Zampini } 9325f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 933b77ba244SStefano Zampini Btype = Bnames ? (size > 1 ? Bnames->mname : Bnames->sname) : Btype; 934b77ba244SStefano Zampini Ctype = Cnames ? (size > 1 ? Cnames->mname : Cnames->sname) : Ctype; 9355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(composedname,sizeof(composedname),"MatProductSetFromOptions_%s_%s_C",((PetscObject)A)->type_name,Btype)); 9365f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetMatProductOperation_Private(A,ptype,symbolic,numeric,destroy,composedname,Ctype)); 9373a40ed3dSBarry Smith PetscFunctionReturn(0); 938e51e0e81SBarry Smith } 939e51e0e81SBarry Smith 9407fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str) 9417fabda3fSAlex Fikl { 94228f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data; 943cb8c8a77SBarry Smith PetscBool matflg; 944b77ba244SStefano Zampini MatShellMatFunctionList matmatA; 9457fabda3fSAlex Fikl 9467fabda3fSAlex Fikl PetscFunctionBegin; 9475f80ce2aSJacob Faibussowitsch CHKERRQ(MatIsShell(B,&matflg)); 948*28b400f6SJacob Faibussowitsch PetscCheck(matflg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix %s not derived from MATSHELL",((PetscObject)B)->type_name); 9497fabda3fSAlex Fikl 9505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps))); 9515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps))); 9527fabda3fSAlex Fikl 953cb8c8a77SBarry Smith if (shellA->ops->copy) { 9545f80ce2aSJacob Faibussowitsch CHKERRQ((*shellA->ops->copy)(A,B,str)); 955cb8c8a77SBarry Smith } 9567fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 9577fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 9580c0fd78eSBarry Smith if (shellA->dshift) { 9590c0fd78eSBarry Smith if (!shellB->dshift) { 9605f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(shellA->dshift,&shellB->dshift)); 9617fabda3fSAlex Fikl } 9625f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(shellA->dshift,shellB->dshift)); 9637fabda3fSAlex Fikl } else { 9645f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shellB->dshift)); 9657fabda3fSAlex Fikl } 9660c0fd78eSBarry Smith if (shellA->left) { 9670c0fd78eSBarry Smith if (!shellB->left) { 9685f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(shellA->left,&shellB->left)); 9697fabda3fSAlex Fikl } 9705f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(shellA->left,shellB->left)); 9717fabda3fSAlex Fikl } else { 9725f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shellB->left)); 9737fabda3fSAlex Fikl } 9740c0fd78eSBarry Smith if (shellA->right) { 9750c0fd78eSBarry Smith if (!shellB->right) { 9765f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(shellA->right,&shellB->right)); 9777fabda3fSAlex Fikl } 9785f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(shellA->right,shellB->right)); 9797fabda3fSAlex Fikl } else { 9805f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shellB->right)); 9817fabda3fSAlex Fikl } 9825f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&shellB->axpy)); 983ef5c7bd2SStefano Zampini shellB->axpy_vscale = 0.0; 984ef5c7bd2SStefano Zampini shellB->axpy_state = 0; 98540e381c3SBarry Smith if (shellA->axpy) { 9865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)shellA->axpy)); 98740e381c3SBarry Smith shellB->axpy = shellA->axpy; 98840e381c3SBarry Smith shellB->axpy_vscale = shellA->axpy_vscale; 989ef5c7bd2SStefano Zampini shellB->axpy_state = shellA->axpy_state; 99040e381c3SBarry Smith } 99145960306SStefano Zampini if (shellA->zrows) { 9925f80ce2aSJacob Faibussowitsch CHKERRQ(ISDuplicate(shellA->zrows,&shellB->zrows)); 99345960306SStefano Zampini if (shellA->zcols) { 9945f80ce2aSJacob Faibussowitsch CHKERRQ(ISDuplicate(shellA->zcols,&shellB->zcols)); 99545960306SStefano Zampini } 9965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(shellA->zvals,&shellB->zvals)); 9975f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(shellA->zvals,shellB->zvals)); 9985f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(shellA->zvals_w,&shellB->zvals_w)); 9995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)shellA->zvals_sct_r)); 10005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)shellA->zvals_sct_c)); 100145960306SStefano Zampini shellB->zvals_sct_r = shellA->zvals_sct_r; 100245960306SStefano Zampini shellB->zvals_sct_c = shellA->zvals_sct_c; 100345960306SStefano Zampini } 1004b77ba244SStefano Zampini 1005b77ba244SStefano Zampini matmatA = shellA->matmat; 1006b77ba244SStefano Zampini if (matmatA) { 1007b77ba244SStefano Zampini while (matmatA->next) { 10085f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetMatProductOperation_Private(B,matmatA->ptype,matmatA->symbolic,matmatA->numeric,matmatA->destroy,matmatA->composedname,matmatA->resultname)); 1009b77ba244SStefano Zampini matmatA = matmatA->next; 1010b77ba244SStefano Zampini } 1011b77ba244SStefano Zampini } 10127fabda3fSAlex Fikl PetscFunctionReturn(0); 10137fabda3fSAlex Fikl } 10147fabda3fSAlex Fikl 1015cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M) 1016cb8c8a77SBarry Smith { 1017cb8c8a77SBarry Smith void *ctx; 1018cb8c8a77SBarry Smith 1019cb8c8a77SBarry Smith PetscFunctionBegin; 10205f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(mat,&ctx)); 10215f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M)); 10225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)(*M),((PetscObject)mat)->type_name)); 1023a4b1380bSStefano Zampini if (op != MAT_DO_NOT_COPY_VALUES) { 10245f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(mat,*M,SAME_NONZERO_PATTERN)); 1025a4b1380bSStefano Zampini } 1026cb8c8a77SBarry Smith PetscFunctionReturn(0); 1027cb8c8a77SBarry Smith } 1028cb8c8a77SBarry Smith 1029dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 1030ef66eb69SBarry Smith { 1031ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 103225578ef6SJed Brown Vec xx; 1033e3f487b0SBarry Smith PetscObjectState instate,outstate; 1034ef66eb69SBarry Smith 1035ef66eb69SBarry Smith PetscFunctionBegin; 1036*28b400f6SJacob Faibussowitsch PetscCheck(shell->ops->mult,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL"); 10375f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellPreZeroRight(A,x,&xx)); 10385f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellPreScaleRight(A,xx,&xx)); 10395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)y, &instate)); 10405f80ce2aSJacob Faibussowitsch CHKERRQ((*shell->ops->mult)(A,xx,y)); 10415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)y, &outstate)); 1042e3f487b0SBarry Smith if (instate == outstate) { 1043e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 10445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateIncrease((PetscObject)y)); 1045e3f487b0SBarry Smith } 10465f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellShiftAndScale(A,xx,y)); 10475f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellPostScaleLeft(A,y)); 10485f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellPostZeroLeft(A,y)); 10499f137db4SBarry Smith 10509f137db4SBarry Smith if (shell->axpy) { 1051ef5c7bd2SStefano Zampini Mat X; 1052ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 1053ef5c7bd2SStefano Zampini 10545f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(shell->axpy,&X)); 10555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)X,&axpy_state)); 10562c71b3e2SJacob Faibussowitsch PetscCheckFalse(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,...)"); 1057b77ba244SStefano Zampini 10585f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(shell->axpy,shell->axpy_right ? NULL : &shell->axpy_right,shell->axpy_left ? NULL : &shell->axpy_left)); 10595f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,shell->axpy_right)); 10605f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(shell->axpy,shell->axpy_right,shell->axpy_left)); 10615f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,shell->axpy_vscale,shell->axpy_left)); 10629f137db4SBarry Smith } 1063ef66eb69SBarry Smith PetscFunctionReturn(0); 1064ef66eb69SBarry Smith } 1065ef66eb69SBarry Smith 10665edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 10675edf6516SJed Brown { 10685edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 10695edf6516SJed Brown 10705edf6516SJed Brown PetscFunctionBegin; 10715edf6516SJed Brown if (y == z) { 10725f80ce2aSJacob Faibussowitsch if (!shell->right_add_work) CHKERRQ(VecDuplicate(z,&shell->right_add_work)); 10735f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,shell->right_add_work)); 10745f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(z,1.0,shell->right_add_work)); 10755edf6516SJed Brown } else { 10765f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,z)); 10775f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(z,1.0,y)); 10785edf6516SJed Brown } 10795edf6516SJed Brown PetscFunctionReturn(0); 10805edf6516SJed Brown } 10815edf6516SJed Brown 108218be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 108318be62a5SSatish Balay { 108418be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 108525578ef6SJed Brown Vec xx; 1086e3f487b0SBarry Smith PetscObjectState instate,outstate; 108718be62a5SSatish Balay 108818be62a5SSatish Balay PetscFunctionBegin; 1089*28b400f6SJacob Faibussowitsch PetscCheck(shell->ops->multtranspose,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL"); 10905f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellPreZeroLeft(A,x,&xx)); 10915f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellPreScaleLeft(A,xx,&xx)); 10925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)y, &instate)); 10935f80ce2aSJacob Faibussowitsch CHKERRQ((*shell->ops->multtranspose)(A,xx,y)); 10945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)y, &outstate)); 1095e3f487b0SBarry Smith if (instate == outstate) { 1096e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 10975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateIncrease((PetscObject)y)); 1098e3f487b0SBarry Smith } 10995f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellShiftAndScale(A,xx,y)); 11005f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellPostScaleRight(A,y)); 11015f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellPostZeroRight(A,y)); 1102350ec83bSStefano Zampini 1103350ec83bSStefano Zampini if (shell->axpy) { 1104ef5c7bd2SStefano Zampini Mat X; 1105ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 1106ef5c7bd2SStefano Zampini 11075f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(shell->axpy,&X)); 11085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)X,&axpy_state)); 11092c71b3e2SJacob Faibussowitsch PetscCheckFalse(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,...)"); 11105f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(shell->axpy,shell->axpy_right ? NULL : &shell->axpy_right,shell->axpy_left ? NULL : &shell->axpy_left)); 11115f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,shell->axpy_left)); 11125f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(shell->axpy,shell->axpy_left,shell->axpy_right)); 11135f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,shell->axpy_vscale,shell->axpy_right)); 1114350ec83bSStefano Zampini } 111518be62a5SSatish Balay PetscFunctionReturn(0); 111618be62a5SSatish Balay } 111718be62a5SSatish Balay 11185edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 11195edf6516SJed Brown { 11205edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 11215edf6516SJed Brown 11225edf6516SJed Brown PetscFunctionBegin; 11235edf6516SJed Brown if (y == z) { 11245f80ce2aSJacob Faibussowitsch if (!shell->left_add_work) CHKERRQ(VecDuplicate(z,&shell->left_add_work)); 11255f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,x,shell->left_add_work)); 11265f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(z,1.0,shell->left_add_work)); 11275edf6516SJed Brown } else { 11285f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,x,z)); 11295f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(z,1.0,y)); 11305edf6516SJed Brown } 11315edf6516SJed Brown PetscFunctionReturn(0); 11325edf6516SJed Brown } 11335edf6516SJed Brown 11340c0fd78eSBarry Smith /* 11350c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 11360c0fd78eSBarry Smith */ 113718be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 113818be62a5SSatish Balay { 113918be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 114018be62a5SSatish Balay 114118be62a5SSatish Balay PetscFunctionBegin; 11420c0fd78eSBarry Smith if (shell->ops->getdiagonal) { 11435f80ce2aSJacob Faibussowitsch CHKERRQ((*shell->ops->getdiagonal)(A,v)); 1144305211d5SBarry 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,...)"); 11455f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(v,shell->vscale)); 11468fe8eb27SJed Brown if (shell->dshift) { 11475f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(v,1.0,shell->dshift)); 114818be62a5SSatish Balay } 11495f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(v,shell->vshift)); 11505f80ce2aSJacob Faibussowitsch if (shell->left) CHKERRQ(VecPointwiseMult(v,v,shell->left)); 11515f80ce2aSJacob Faibussowitsch if (shell->right) CHKERRQ(VecPointwiseMult(v,v,shell->right)); 115245960306SStefano Zampini if (shell->zrows) { 11535f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE)); 11545f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE)); 115545960306SStefano Zampini } 115681c02519SBarry Smith if (shell->axpy) { 1157ef5c7bd2SStefano Zampini Mat X; 1158ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 1159ef5c7bd2SStefano Zampini 11605f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(shell->axpy,&X)); 11615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)X,&axpy_state)); 11622c71b3e2SJacob Faibussowitsch PetscCheckFalse(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,...)"); 11635f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(shell->axpy,NULL,shell->axpy_left ? NULL : &shell->axpy_left)); 11645f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(shell->axpy,shell->axpy_left)); 11655f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(v,shell->axpy_vscale,shell->axpy_left)); 116681c02519SBarry Smith } 116718be62a5SSatish Balay PetscFunctionReturn(0); 116818be62a5SSatish Balay } 116918be62a5SSatish Balay 1170f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 1171ef66eb69SBarry Smith { 1172ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 1173789d8953SBarry Smith PetscBool flg; 1174b24ad042SBarry Smith 1175ef66eb69SBarry Smith PetscFunctionBegin; 11765f80ce2aSJacob Faibussowitsch CHKERRQ(MatHasCongruentLayouts(Y,&flg)); 1177*28b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot shift shell matrix if it is not congruent"); 11780c0fd78eSBarry Smith if (shell->left || shell->right) { 11798fe8eb27SJed Brown if (!shell->dshift) { 11805f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift)); 11815f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(shell->dshift,a)); 11820c0fd78eSBarry Smith } else { 11835f80ce2aSJacob Faibussowitsch if (shell->left) CHKERRQ(VecPointwiseMult(shell->dshift,shell->dshift,shell->left)); 11845f80ce2aSJacob Faibussowitsch if (shell->right) CHKERRQ(VecPointwiseMult(shell->dshift,shell->dshift,shell->right)); 11855f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(shell->dshift,a)); 11860c0fd78eSBarry Smith } 11875f80ce2aSJacob Faibussowitsch if (shell->left) CHKERRQ(VecPointwiseDivide(shell->dshift,shell->dshift,shell->left)); 11885f80ce2aSJacob Faibussowitsch if (shell->right) CHKERRQ(VecPointwiseDivide(shell->dshift,shell->dshift,shell->right)); 11898fe8eb27SJed Brown } else shell->vshift += a; 119045960306SStefano Zampini if (shell->zrows) { 11915f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(shell->zvals,a)); 119245960306SStefano Zampini } 1193ef66eb69SBarry Smith PetscFunctionReturn(0); 1194ef66eb69SBarry Smith } 1195ef66eb69SBarry Smith 1196b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s) 11976464bf51SAlex Fikl { 11986464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 11996464bf51SAlex Fikl 12006464bf51SAlex Fikl PetscFunctionBegin; 12015f80ce2aSJacob Faibussowitsch if (!shell->dshift) CHKERRQ(VecDuplicate(D,&shell->dshift)); 12020c0fd78eSBarry Smith if (shell->left || shell->right) { 12035f80ce2aSJacob Faibussowitsch if (!shell->right_work) CHKERRQ(VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work)); 12049f137db4SBarry Smith if (shell->left && shell->right) { 12055f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseDivide(shell->right_work,D,shell->left)); 12065f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseDivide(shell->right_work,shell->right_work,shell->right)); 12079f137db4SBarry Smith } else if (shell->left) { 12085f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseDivide(shell->right_work,D,shell->left)); 12099f137db4SBarry Smith } else { 12105f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseDivide(shell->right_work,D,shell->right)); 12119f137db4SBarry Smith } 12125f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(shell->dshift,s,shell->right_work)); 12130c0fd78eSBarry Smith } else { 12145f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(shell->dshift,s,D)); 1215b253ae0bS“Marek } 1216b253ae0bS“Marek PetscFunctionReturn(0); 1217b253ae0bS“Marek } 1218b253ae0bS“Marek 1219b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 1220b253ae0bS“Marek { 122145960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 1222b253ae0bS“Marek Vec d; 1223789d8953SBarry Smith PetscBool flg; 1224b253ae0bS“Marek 1225b253ae0bS“Marek PetscFunctionBegin; 12265f80ce2aSJacob Faibussowitsch CHKERRQ(MatHasCongruentLayouts(A,&flg)); 1227*28b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot diagonal set or shift shell matrix if it is not congruent"); 1228b253ae0bS“Marek if (ins == INSERT_VALUES) { 1229*28b400f6SJacob Faibussowitsch PetscCheck(A->ops->getdiagonal,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first"); 12305f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(D,&d)); 12315f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(A,d)); 12325f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalSet_Shell_Private(A,d,-1.)); 12335f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalSet_Shell_Private(A,D,1.)); 12345f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d)); 123545960306SStefano Zampini if (shell->zrows) { 12365f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(D,shell->zvals)); 123745960306SStefano Zampini } 1238b253ae0bS“Marek } else { 12395f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalSet_Shell_Private(A,D,1.)); 124045960306SStefano Zampini if (shell->zrows) { 12415f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(shell->zvals,1.0,D)); 124245960306SStefano Zampini } 12436464bf51SAlex Fikl } 12446464bf51SAlex Fikl PetscFunctionReturn(0); 12456464bf51SAlex Fikl } 12466464bf51SAlex Fikl 1247f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 1248ef66eb69SBarry Smith { 1249ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 1250b24ad042SBarry Smith 1251ef66eb69SBarry Smith PetscFunctionBegin; 1252f4df32b1SMatthew Knepley shell->vscale *= a; 12530c0fd78eSBarry Smith shell->vshift *= a; 12542205254eSKarl Rupp if (shell->dshift) { 12555f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(shell->dshift,a)); 12560c0fd78eSBarry Smith } 125781c02519SBarry Smith shell->axpy_vscale *= a; 125845960306SStefano Zampini if (shell->zrows) { 12595f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(shell->zvals,a)); 126045960306SStefano Zampini } 12618fe8eb27SJed Brown PetscFunctionReturn(0); 126218be62a5SSatish Balay } 12638fe8eb27SJed Brown 12648fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 12658fe8eb27SJed Brown { 12668fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 12678fe8eb27SJed Brown 12688fe8eb27SJed Brown PetscFunctionBegin; 12698fe8eb27SJed Brown if (left) { 12700c0fd78eSBarry Smith if (!shell->left) { 12715f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(left,&shell->left)); 12725f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(left,shell->left)); 12730c0fd78eSBarry Smith } else { 12745f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->left,shell->left,left)); 127518be62a5SSatish Balay } 127645960306SStefano Zampini if (shell->zrows) { 12775f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->zvals,shell->zvals,left)); 127845960306SStefano Zampini } 1279ef66eb69SBarry Smith } 12808fe8eb27SJed Brown if (right) { 12810c0fd78eSBarry Smith if (!shell->right) { 12825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(right,&shell->right)); 12835f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(right,shell->right)); 12840c0fd78eSBarry Smith } else { 12855f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->right,shell->right,right)); 12868fe8eb27SJed Brown } 128745960306SStefano Zampini if (shell->zrows) { 128845960306SStefano Zampini if (!shell->left_work) { 12895f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(Y,NULL,&shell->left_work)); 129045960306SStefano Zampini } 12915f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(shell->zvals_w,1.0)); 12925f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 12935f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 12945f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(shell->zvals,shell->zvals,shell->zvals_w)); 129545960306SStefano Zampini } 12968fe8eb27SJed Brown } 1297ef5c7bd2SStefano Zampini if (shell->axpy) { 12985f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(shell->axpy,left,right)); 1299ef5c7bd2SStefano Zampini } 1300ef66eb69SBarry Smith PetscFunctionReturn(0); 1301ef66eb69SBarry Smith } 1302ef66eb69SBarry Smith 1303dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 1304ef66eb69SBarry Smith { 1305ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 1306ef66eb69SBarry Smith 1307ef66eb69SBarry Smith PetscFunctionBegin; 13088fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 1309ef66eb69SBarry Smith shell->vshift = 0.0; 1310ef66eb69SBarry Smith shell->vscale = 1.0; 1311ef5c7bd2SStefano Zampini shell->axpy_vscale = 0.0; 1312ef5c7bd2SStefano Zampini shell->axpy_state = 0; 13135f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->dshift)); 13145f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->left)); 13155f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->right)); 13165f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&shell->axpy)); 13175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->axpy_left)); 13185f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&shell->axpy_right)); 13195f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&shell->zvals_sct_c)); 13205f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&shell->zvals_sct_r)); 13215f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&shell->zrows)); 13225f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&shell->zcols)); 1323ef66eb69SBarry Smith } 1324ef66eb69SBarry Smith PetscFunctionReturn(0); 1325ef66eb69SBarry Smith } 1326ef66eb69SBarry Smith 13273b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 13283b49f96aSBarry Smith { 13293b49f96aSBarry Smith PetscFunctionBegin; 13303b49f96aSBarry Smith *missing = PETSC_FALSE; 13313b49f96aSBarry Smith PetscFunctionReturn(0); 13323b49f96aSBarry Smith } 13333b49f96aSBarry Smith 13349f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str) 13359f137db4SBarry Smith { 13369f137db4SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 13379f137db4SBarry Smith 13389f137db4SBarry Smith PetscFunctionBegin; 1339ef5c7bd2SStefano Zampini if (X == Y) { 13405f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(Y,1.0 + a)); 1341ef5c7bd2SStefano Zampini PetscFunctionReturn(0); 1342ef5c7bd2SStefano Zampini } 1343ef5c7bd2SStefano Zampini if (!shell->axpy) { 13445f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvertFrom_Shell(X,MATSHELL,MAT_INITIAL_MATRIX,&shell->axpy)); 13459f137db4SBarry Smith shell->axpy_vscale = a; 13465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)X,&shell->axpy_state)); 1347ef5c7bd2SStefano Zampini } else { 13485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(shell->axpy,a/shell->axpy_vscale,X,str)); 1349ef5c7bd2SStefano Zampini } 13509f137db4SBarry Smith PetscFunctionReturn(0); 13519f137db4SBarry Smith } 13529f137db4SBarry Smith 1353f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL, 1354f4259b30SLisandro Dalcin NULL, 1355f4259b30SLisandro Dalcin NULL, 1356f4259b30SLisandro Dalcin NULL, 13570c0fd78eSBarry Smith /* 4*/ MatMultAdd_Shell, 1358f4259b30SLisandro Dalcin NULL, 13590c0fd78eSBarry Smith MatMultTransposeAdd_Shell, 1360f4259b30SLisandro Dalcin NULL, 1361f4259b30SLisandro Dalcin NULL, 1362f4259b30SLisandro Dalcin NULL, 1363f4259b30SLisandro Dalcin /*10*/ NULL, 1364f4259b30SLisandro Dalcin NULL, 1365f4259b30SLisandro Dalcin NULL, 1366f4259b30SLisandro Dalcin NULL, 1367f4259b30SLisandro Dalcin NULL, 1368f4259b30SLisandro Dalcin /*15*/ NULL, 1369f4259b30SLisandro Dalcin NULL, 1370f4259b30SLisandro Dalcin NULL, 13718fe8eb27SJed Brown MatDiagonalScale_Shell, 1372f4259b30SLisandro Dalcin NULL, 1373f4259b30SLisandro Dalcin /*20*/ NULL, 1374ef66eb69SBarry Smith MatAssemblyEnd_Shell, 1375f4259b30SLisandro Dalcin NULL, 1376f4259b30SLisandro Dalcin NULL, 137745960306SStefano Zampini /*24*/ MatZeroRows_Shell, 1378f4259b30SLisandro Dalcin NULL, 1379f4259b30SLisandro Dalcin NULL, 1380f4259b30SLisandro Dalcin NULL, 1381f4259b30SLisandro Dalcin NULL, 1382f4259b30SLisandro Dalcin /*29*/ NULL, 1383f4259b30SLisandro Dalcin NULL, 1384f4259b30SLisandro Dalcin NULL, 1385f4259b30SLisandro Dalcin NULL, 1386f4259b30SLisandro Dalcin NULL, 1387cb8c8a77SBarry Smith /*34*/ MatDuplicate_Shell, 1388f4259b30SLisandro Dalcin NULL, 1389f4259b30SLisandro Dalcin NULL, 1390f4259b30SLisandro Dalcin NULL, 1391f4259b30SLisandro Dalcin NULL, 13929f137db4SBarry Smith /*39*/ MatAXPY_Shell, 1393f4259b30SLisandro Dalcin NULL, 1394f4259b30SLisandro Dalcin NULL, 1395f4259b30SLisandro Dalcin NULL, 1396cb8c8a77SBarry Smith MatCopy_Shell, 1397f4259b30SLisandro Dalcin /*44*/ NULL, 1398ef66eb69SBarry Smith MatScale_Shell, 1399ef66eb69SBarry Smith MatShift_Shell, 14006464bf51SAlex Fikl MatDiagonalSet_Shell, 140145960306SStefano Zampini MatZeroRowsColumns_Shell, 1402f4259b30SLisandro Dalcin /*49*/ NULL, 1403f4259b30SLisandro Dalcin NULL, 1404f4259b30SLisandro Dalcin NULL, 1405f4259b30SLisandro Dalcin NULL, 1406f4259b30SLisandro Dalcin NULL, 1407f4259b30SLisandro Dalcin /*54*/ NULL, 1408f4259b30SLisandro Dalcin NULL, 1409f4259b30SLisandro Dalcin NULL, 1410f4259b30SLisandro Dalcin NULL, 1411f4259b30SLisandro Dalcin NULL, 1412f4259b30SLisandro Dalcin /*59*/ NULL, 1413b9b97703SBarry Smith MatDestroy_Shell, 1414f4259b30SLisandro Dalcin NULL, 1415251fad3fSStefano Zampini MatConvertFrom_Shell, 1416f4259b30SLisandro Dalcin NULL, 1417f4259b30SLisandro Dalcin /*64*/ NULL, 1418f4259b30SLisandro Dalcin NULL, 1419f4259b30SLisandro Dalcin NULL, 1420f4259b30SLisandro Dalcin NULL, 1421f4259b30SLisandro Dalcin NULL, 1422f4259b30SLisandro Dalcin /*69*/ NULL, 1423f4259b30SLisandro Dalcin NULL, 1424c87e5d42SMatthew Knepley MatConvert_Shell, 1425f4259b30SLisandro Dalcin NULL, 1426f4259b30SLisandro Dalcin NULL, 1427f4259b30SLisandro Dalcin /*74*/ NULL, 1428f4259b30SLisandro Dalcin NULL, 1429f4259b30SLisandro Dalcin NULL, 1430f4259b30SLisandro Dalcin NULL, 1431f4259b30SLisandro Dalcin NULL, 1432f4259b30SLisandro Dalcin /*79*/ NULL, 1433f4259b30SLisandro Dalcin NULL, 1434f4259b30SLisandro Dalcin NULL, 1435f4259b30SLisandro Dalcin NULL, 1436f4259b30SLisandro Dalcin NULL, 1437f4259b30SLisandro Dalcin /*84*/ NULL, 1438f4259b30SLisandro Dalcin NULL, 1439f4259b30SLisandro Dalcin NULL, 1440f4259b30SLisandro Dalcin NULL, 1441f4259b30SLisandro Dalcin NULL, 1442f4259b30SLisandro Dalcin /*89*/ NULL, 1443f4259b30SLisandro Dalcin NULL, 1444f4259b30SLisandro Dalcin NULL, 1445f4259b30SLisandro Dalcin NULL, 1446f4259b30SLisandro Dalcin NULL, 1447f4259b30SLisandro Dalcin /*94*/ NULL, 1448f4259b30SLisandro Dalcin NULL, 1449f4259b30SLisandro Dalcin NULL, 1450f4259b30SLisandro Dalcin NULL, 1451f4259b30SLisandro Dalcin NULL, 1452f4259b30SLisandro Dalcin /*99*/ NULL, 1453f4259b30SLisandro Dalcin NULL, 1454f4259b30SLisandro Dalcin NULL, 1455f4259b30SLisandro Dalcin NULL, 1456f4259b30SLisandro Dalcin NULL, 1457f4259b30SLisandro Dalcin /*104*/ NULL, 1458f4259b30SLisandro Dalcin NULL, 1459f4259b30SLisandro Dalcin NULL, 1460f4259b30SLisandro Dalcin NULL, 1461f4259b30SLisandro Dalcin NULL, 1462f4259b30SLisandro Dalcin /*109*/ NULL, 1463f4259b30SLisandro Dalcin NULL, 1464f4259b30SLisandro Dalcin NULL, 1465f4259b30SLisandro Dalcin NULL, 14663b49f96aSBarry Smith MatMissingDiagonal_Shell, 1467f4259b30SLisandro Dalcin /*114*/ NULL, 1468f4259b30SLisandro Dalcin NULL, 1469f4259b30SLisandro Dalcin NULL, 1470f4259b30SLisandro Dalcin NULL, 1471f4259b30SLisandro Dalcin NULL, 1472f4259b30SLisandro Dalcin /*119*/ NULL, 1473f4259b30SLisandro Dalcin NULL, 1474f4259b30SLisandro Dalcin NULL, 1475f4259b30SLisandro Dalcin NULL, 1476f4259b30SLisandro Dalcin NULL, 1477f4259b30SLisandro Dalcin /*124*/ NULL, 1478f4259b30SLisandro Dalcin NULL, 1479f4259b30SLisandro Dalcin NULL, 1480f4259b30SLisandro Dalcin NULL, 1481f4259b30SLisandro Dalcin NULL, 1482f4259b30SLisandro Dalcin /*129*/ NULL, 1483f4259b30SLisandro Dalcin NULL, 1484f4259b30SLisandro Dalcin NULL, 1485f4259b30SLisandro Dalcin NULL, 1486f4259b30SLisandro Dalcin NULL, 1487f4259b30SLisandro Dalcin /*134*/ NULL, 1488f4259b30SLisandro Dalcin NULL, 1489f4259b30SLisandro Dalcin NULL, 1490f4259b30SLisandro Dalcin NULL, 1491f4259b30SLisandro Dalcin NULL, 1492f4259b30SLisandro Dalcin /*139*/ NULL, 1493f4259b30SLisandro Dalcin NULL, 1494d70f29a3SPierre Jolivet NULL, 1495d70f29a3SPierre Jolivet NULL, 1496d70f29a3SPierre Jolivet NULL, 1497d70f29a3SPierre Jolivet /*144*/ NULL, 1498d70f29a3SPierre Jolivet NULL, 1499d70f29a3SPierre Jolivet NULL, 1500f4259b30SLisandro Dalcin NULL 15013964eb88SJed Brown }; 1502273d9f13SBarry Smith 1503789d8953SBarry Smith PetscErrorCode MatShellSetContext_Shell(Mat mat,void *ctx) 1504789d8953SBarry Smith { 1505789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1506789d8953SBarry Smith 1507789d8953SBarry Smith PetscFunctionBegin; 1508789d8953SBarry Smith shell->ctx = ctx; 1509789d8953SBarry Smith PetscFunctionReturn(0); 1510789d8953SBarry Smith } 1511789d8953SBarry Smith 1512db77db73SJeremy L Thompson static PetscErrorCode MatShellSetVecType_Shell(Mat mat,VecType vtype) 1513db77db73SJeremy L Thompson { 1514db77db73SJeremy L Thompson PetscFunctionBegin; 15155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(mat->defaultvectype)); 15165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(vtype,(char**)&mat->defaultvectype)); 1517db77db73SJeremy L Thompson PetscFunctionReturn(0); 1518db77db73SJeremy L Thompson } 1519db77db73SJeremy L Thompson 1520789d8953SBarry Smith PetscErrorCode MatShellSetManageScalingShifts_Shell(Mat A) 1521789d8953SBarry Smith { 1522789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 1523789d8953SBarry Smith 1524789d8953SBarry Smith PetscFunctionBegin; 1525789d8953SBarry Smith shell->managescalingshifts = PETSC_FALSE; 1526789d8953SBarry Smith A->ops->diagonalset = NULL; 1527789d8953SBarry Smith A->ops->diagonalscale = NULL; 1528789d8953SBarry Smith A->ops->scale = NULL; 1529789d8953SBarry Smith A->ops->shift = NULL; 1530789d8953SBarry Smith A->ops->axpy = NULL; 1531789d8953SBarry Smith PetscFunctionReturn(0); 1532789d8953SBarry Smith } 1533789d8953SBarry Smith 1534789d8953SBarry Smith PetscErrorCode MatShellSetOperation_Shell(Mat mat,MatOperation op,void (*f)(void)) 1535789d8953SBarry Smith { 1536feb237baSPierre Jolivet Mat_Shell *shell = (Mat_Shell*)mat->data; 1537789d8953SBarry Smith 1538789d8953SBarry Smith PetscFunctionBegin; 1539789d8953SBarry Smith switch (op) { 1540789d8953SBarry Smith case MATOP_DESTROY: 1541789d8953SBarry Smith shell->ops->destroy = (PetscErrorCode (*)(Mat))f; 1542789d8953SBarry Smith break; 1543789d8953SBarry Smith case MATOP_VIEW: 1544789d8953SBarry Smith if (!mat->ops->viewnative) { 1545789d8953SBarry Smith mat->ops->viewnative = mat->ops->view; 1546789d8953SBarry Smith } 1547789d8953SBarry Smith mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 1548789d8953SBarry Smith break; 1549789d8953SBarry Smith case MATOP_COPY: 1550789d8953SBarry Smith shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f; 1551789d8953SBarry Smith break; 1552789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1553789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1554789d8953SBarry Smith case MATOP_SHIFT: 1555789d8953SBarry Smith case MATOP_SCALE: 1556789d8953SBarry Smith case MATOP_AXPY: 1557789d8953SBarry Smith case MATOP_ZERO_ROWS: 1558789d8953SBarry Smith case MATOP_ZERO_ROWS_COLUMNS: 1559*28b400f6SJacob Faibussowitsch PetscCheck(!shell->managescalingshifts,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()"); 1560789d8953SBarry Smith (((void(**)(void))mat->ops)[op]) = f; 1561789d8953SBarry Smith break; 1562789d8953SBarry Smith case MATOP_GET_DIAGONAL: 1563789d8953SBarry Smith if (shell->managescalingshifts) { 1564789d8953SBarry Smith shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1565789d8953SBarry Smith mat->ops->getdiagonal = MatGetDiagonal_Shell; 1566789d8953SBarry Smith } else { 1567789d8953SBarry Smith shell->ops->getdiagonal = NULL; 1568789d8953SBarry Smith mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1569789d8953SBarry Smith } 1570789d8953SBarry Smith break; 1571789d8953SBarry Smith case MATOP_MULT: 1572789d8953SBarry Smith if (shell->managescalingshifts) { 1573789d8953SBarry Smith shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1574789d8953SBarry Smith mat->ops->mult = MatMult_Shell; 1575789d8953SBarry Smith } else { 1576789d8953SBarry Smith shell->ops->mult = NULL; 1577789d8953SBarry Smith mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1578789d8953SBarry Smith } 1579789d8953SBarry Smith break; 1580789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 1581789d8953SBarry Smith if (shell->managescalingshifts) { 1582789d8953SBarry Smith shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1583789d8953SBarry Smith mat->ops->multtranspose = MatMultTranspose_Shell; 1584789d8953SBarry Smith } else { 1585789d8953SBarry Smith shell->ops->multtranspose = NULL; 1586789d8953SBarry Smith mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1587789d8953SBarry Smith } 1588789d8953SBarry Smith break; 1589789d8953SBarry Smith default: 1590789d8953SBarry Smith (((void(**)(void))mat->ops)[op]) = f; 1591789d8953SBarry Smith break; 1592789d8953SBarry Smith } 1593789d8953SBarry Smith PetscFunctionReturn(0); 1594789d8953SBarry Smith } 1595789d8953SBarry Smith 1596789d8953SBarry Smith PetscErrorCode MatShellGetOperation_Shell(Mat mat,MatOperation op,void(**f)(void)) 1597789d8953SBarry Smith { 1598789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1599789d8953SBarry Smith 1600789d8953SBarry Smith PetscFunctionBegin; 1601789d8953SBarry Smith switch (op) { 1602789d8953SBarry Smith case MATOP_DESTROY: 1603789d8953SBarry Smith *f = (void (*)(void))shell->ops->destroy; 1604789d8953SBarry Smith break; 1605789d8953SBarry Smith case MATOP_VIEW: 1606789d8953SBarry Smith *f = (void (*)(void))mat->ops->view; 1607789d8953SBarry Smith break; 1608789d8953SBarry Smith case MATOP_COPY: 1609789d8953SBarry Smith *f = (void (*)(void))shell->ops->copy; 1610789d8953SBarry Smith break; 1611789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1612789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1613789d8953SBarry Smith case MATOP_SHIFT: 1614789d8953SBarry Smith case MATOP_SCALE: 1615789d8953SBarry Smith case MATOP_AXPY: 1616789d8953SBarry Smith case MATOP_ZERO_ROWS: 1617789d8953SBarry Smith case MATOP_ZERO_ROWS_COLUMNS: 1618789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1619789d8953SBarry Smith break; 1620789d8953SBarry Smith case MATOP_GET_DIAGONAL: 1621789d8953SBarry Smith if (shell->ops->getdiagonal) 1622789d8953SBarry Smith *f = (void (*)(void))shell->ops->getdiagonal; 1623789d8953SBarry Smith else 1624789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1625789d8953SBarry Smith break; 1626789d8953SBarry Smith case MATOP_MULT: 1627789d8953SBarry Smith if (shell->ops->mult) 1628789d8953SBarry Smith *f = (void (*)(void))shell->ops->mult; 1629789d8953SBarry Smith else 1630789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1631789d8953SBarry Smith break; 1632789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 1633789d8953SBarry Smith if (shell->ops->multtranspose) 1634789d8953SBarry Smith *f = (void (*)(void))shell->ops->multtranspose; 1635789d8953SBarry Smith else 1636789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1637789d8953SBarry Smith break; 1638789d8953SBarry Smith default: 1639789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1640789d8953SBarry Smith } 1641789d8953SBarry Smith PetscFunctionReturn(0); 1642789d8953SBarry Smith } 1643789d8953SBarry Smith 16440bad9183SKris Buschelman /*MC 1645fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 16460bad9183SKris Buschelman 16470bad9183SKris Buschelman Level: advanced 16480bad9183SKris Buschelman 16490c0fd78eSBarry Smith .seealso: MatCreateShell() 16500bad9183SKris Buschelman M*/ 16510bad9183SKris Buschelman 16528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 1653273d9f13SBarry Smith { 1654273d9f13SBarry Smith Mat_Shell *b; 1655273d9f13SBarry Smith 1656273d9f13SBarry Smith PetscFunctionBegin; 16575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps))); 1658273d9f13SBarry Smith 16595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(A,&b)); 1660273d9f13SBarry Smith A->data = (void*)b; 1661273d9f13SBarry Smith 1662f4259b30SLisandro Dalcin b->ctx = NULL; 1663ef66eb69SBarry Smith b->vshift = 0.0; 1664ef66eb69SBarry Smith b->vscale = 1.0; 16650c0fd78eSBarry Smith b->managescalingshifts = PETSC_TRUE; 1666273d9f13SBarry Smith A->assembled = PETSC_TRUE; 1667210f0121SBarry Smith A->preallocated = PETSC_FALSE; 16682205254eSKarl Rupp 16695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatShellGetContext_C",MatShellGetContext_Shell)); 16705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatShellSetContext_C",MatShellSetContext_Shell)); 16715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatShellSetVecType_C",MatShellSetVecType_Shell)); 16725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatShellSetManageScalingShifts_C",MatShellSetManageScalingShifts_Shell)); 16735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatShellSetOperation_C",MatShellSetOperation_Shell)); 16745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatShellGetOperation_C",MatShellGetOperation_Shell)); 16755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatShellSetMatProductOperation_C",MatShellSetMatProductOperation_Shell)); 16765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)A,MATSHELL)); 1677273d9f13SBarry Smith PetscFunctionReturn(0); 1678273d9f13SBarry Smith } 1679e51e0e81SBarry Smith 16804b828684SBarry Smith /*@C 1681052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 1682ff756334SLois Curfman McInnes private data storage format. 1683e51e0e81SBarry Smith 1684d083f849SBarry Smith Collective 1685c7fcc2eaSBarry Smith 1686e51e0e81SBarry Smith Input Parameters: 1687c7fcc2eaSBarry Smith + comm - MPI communicator 1688c7fcc2eaSBarry Smith . m - number of local rows (must be given) 1689c7fcc2eaSBarry Smith . n - number of local columns (must be given) 1690c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 1691c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 1692c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 1693e51e0e81SBarry Smith 1694ff756334SLois Curfman McInnes Output Parameter: 169544cd7ae7SLois Curfman McInnes . A - the matrix 1696e51e0e81SBarry Smith 1697ff2fd236SBarry Smith Level: advanced 1698ff2fd236SBarry Smith 1699f39d1f56SLois Curfman McInnes Usage: 17005bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1701f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 1702c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1703f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 1704f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 1705f39d1f56SLois Curfman McInnes 1706ff756334SLois Curfman McInnes Notes: 1707ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 1708ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 1709ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 1710e51e0e81SBarry Smith 171195452b02SPatrick Sanan Fortran Notes: 171295452b02SPatrick Sanan To use this from Fortran with a ctx you must write an interface definition for this 1713daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 1714daf670e6SBarry Smith in as the ctx argument. 1715f38a8d56SBarry Smith 1716f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 1717f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 1718645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 1719645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 1720645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 1721645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 1722645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 1723645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 1724645985a0SLois Curfman McInnes For example, 1725f39d1f56SLois Curfman McInnes 1726f39d1f56SLois Curfman McInnes $ 1727f39d1f56SLois Curfman McInnes $ Vec x, y 17285bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1729645985a0SLois Curfman McInnes $ Mat A 1730f39d1f56SLois Curfman McInnes $ 1731c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 1732c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 1733f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 1734c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 1735c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 1736c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1737645985a0SLois Curfman McInnes $ MatMult(A,x,y); 173845960306SStefano Zampini $ MatDestroy(&A); 173945960306SStefano Zampini $ VecDestroy(&y); 174045960306SStefano Zampini $ VecDestroy(&x); 1741645985a0SLois Curfman McInnes $ 1742e51e0e81SBarry Smith 174345960306SStefano Zampini MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), MatScale(), MatZeroRows() and MatZeroRowsColumns() internally, so these 17449b53d762SBarry Smith operations cannot be overwritten unless MatShellSetManageScalingShifts() is called. 17459b53d762SBarry Smith 17460c0fd78eSBarry Smith For rectangular matrices do all the scalings and shifts make sense? 17470c0fd78eSBarry Smith 174895452b02SPatrick Sanan Developers Notes: 174995452b02SPatrick Sanan Regarding shifting and scaling. The general form is 17500c0fd78eSBarry Smith 17510c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 17520c0fd78eSBarry Smith 17530c0fd78eSBarry Smith The order you apply the operations is important. For example if you have a dshift then 17540c0fd78eSBarry Smith apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift 17550c0fd78eSBarry Smith you get s*vscale*A + diag(shift) 17560c0fd78eSBarry Smith 17570c0fd78eSBarry Smith A is the user provided function. 17580c0fd78eSBarry Smith 1759ad2f5c3fSBarry Smith KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for 1760ad2f5c3fSBarry Smith for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger 1761ad2f5c3fSBarry Smith an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat); 1762ad2f5c3fSBarry Smith each time the MATSHELL matrix has changed. 1763ad2f5c3fSBarry Smith 17647301b172SPierre Jolivet Matrix product operations (i.e. MatMat, MatTransposeMat etc) can be specified using MatShellSetMatProductOperation() 1765b77ba244SStefano Zampini 1766ad2f5c3fSBarry Smith Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided 1767ad2f5c3fSBarry Smith with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale(). 1768ad2f5c3fSBarry Smith 1769b77ba244SStefano Zampini .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts(), MatShellSetMatProductOperation() 1770e51e0e81SBarry Smith @*/ 17717087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 1772e51e0e81SBarry Smith { 17733a40ed3dSBarry Smith PetscFunctionBegin; 17745f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,A)); 17755f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(*A,m,n,M,N)); 17765f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(*A,MATSHELL)); 17775f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetContext(*A,ctx)); 17785f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(*A)); 1779273d9f13SBarry Smith PetscFunctionReturn(0); 1780c7fcc2eaSBarry Smith } 1781c7fcc2eaSBarry Smith 1782c6866cfdSSatish Balay /*@ 1783273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 1784c7fcc2eaSBarry Smith 17853f9fe445SBarry Smith Logically Collective on Mat 1786c7fcc2eaSBarry Smith 1787273d9f13SBarry Smith Input Parameters: 1788273d9f13SBarry Smith + mat - the shell matrix 1789273d9f13SBarry Smith - ctx - the context 1790273d9f13SBarry Smith 1791273d9f13SBarry Smith Level: advanced 1792273d9f13SBarry Smith 179395452b02SPatrick Sanan Fortran Notes: 179495452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 1795daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 1796273d9f13SBarry Smith 1797273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 17980bc0a719SBarry Smith @*/ 17997087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 1800273d9f13SBarry Smith { 1801273d9f13SBarry Smith PetscFunctionBegin; 18020700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 18035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(mat,"MatShellSetContext_C",(Mat,void*),(mat,ctx))); 18043a40ed3dSBarry Smith PetscFunctionReturn(0); 1805e51e0e81SBarry Smith } 1806e51e0e81SBarry Smith 1807db77db73SJeremy L Thompson /*@C 1808db77db73SJeremy L Thompson MatShellSetVecType - Sets the type of Vec returned by MatCreateVecs() 1809db77db73SJeremy L Thompson 1810db77db73SJeremy L Thompson Logically collective 1811db77db73SJeremy L Thompson 1812db77db73SJeremy L Thompson Input Parameters: 1813db77db73SJeremy L Thompson + mat - the shell matrix 1814db77db73SJeremy L Thompson - vtype - type to use for creating vectors 1815db77db73SJeremy L Thompson 1816db77db73SJeremy L Thompson Notes: 1817db77db73SJeremy L Thompson 1818db77db73SJeremy L Thompson Level: advanced 1819db77db73SJeremy L Thompson 1820db77db73SJeremy L Thompson .seealso: MatCreateVecs() 1821db77db73SJeremy L Thompson @*/ 1822db77db73SJeremy L Thompson PetscErrorCode MatShellSetVecType(Mat mat,VecType vtype) 1823db77db73SJeremy L Thompson { 1824db77db73SJeremy L Thompson PetscFunctionBegin; 18255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(mat,"MatShellSetVecType_C",(Mat,VecType),(mat,vtype))); 1826db77db73SJeremy L Thompson PetscFunctionReturn(0); 1827db77db73SJeremy L Thompson } 1828db77db73SJeremy L Thompson 18290c0fd78eSBarry Smith /*@ 18300c0fd78eSBarry Smith MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately 18310c0fd78eSBarry Smith after MatCreateShell() 18320c0fd78eSBarry Smith 18330c0fd78eSBarry Smith Logically Collective on Mat 18340c0fd78eSBarry Smith 18350c0fd78eSBarry Smith Input Parameter: 18360c0fd78eSBarry Smith . mat - the shell matrix 18370c0fd78eSBarry Smith 18380c0fd78eSBarry Smith Level: advanced 18390c0fd78eSBarry Smith 18400c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation() 18410c0fd78eSBarry Smith @*/ 18420c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A) 18430c0fd78eSBarry Smith { 18440c0fd78eSBarry Smith PetscFunctionBegin; 18450c0fd78eSBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 18465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(A,"MatShellSetManageScalingShifts_C",(Mat),(A))); 18470c0fd78eSBarry Smith PetscFunctionReturn(0); 18480c0fd78eSBarry Smith } 18490c0fd78eSBarry Smith 1850c16cb8f2SBarry Smith /*@C 1851f3b1f45cSBarry Smith MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function. 1852f3b1f45cSBarry Smith 1853f3b1f45cSBarry Smith Logically Collective on Mat 1854f3b1f45cSBarry Smith 1855f3b1f45cSBarry Smith Input Parameters: 1856f3b1f45cSBarry Smith + mat - the shell matrix 1857f3b1f45cSBarry Smith . f - the function 1858f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1859f3b1f45cSBarry Smith - ctx - an optional context for the function 1860f3b1f45cSBarry Smith 1861f3b1f45cSBarry Smith Output Parameter: 1862f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1863f3b1f45cSBarry Smith 1864f3b1f45cSBarry Smith Options Database: 1865f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1866f3b1f45cSBarry Smith 1867f3b1f45cSBarry Smith Level: advanced 1868f3b1f45cSBarry Smith 186995452b02SPatrick Sanan Fortran Notes: 187095452b02SPatrick Sanan Not supported from Fortran 1871f3b1f45cSBarry Smith 1872f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose() 1873f3b1f45cSBarry Smith @*/ 1874f3b1f45cSBarry Smith PetscErrorCode MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1875f3b1f45cSBarry Smith { 1876f3b1f45cSBarry Smith PetscInt m,n; 1877f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1878f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 187974e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1880f3b1f45cSBarry Smith 1881f3b1f45cSBarry Smith PetscFunctionBegin; 1882f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 18835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v)); 18845f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(mat,&m,&n)); 18855f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf)); 18865f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetFunction(mf,f,ctx)); 18875f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetBase(mf,base,NULL)); 1888f3b1f45cSBarry Smith 18895f80ce2aSJacob Faibussowitsch CHKERRQ(MatComputeOperator(mf,MATAIJ,&Dmf)); 18905f80ce2aSJacob Faibussowitsch CHKERRQ(MatComputeOperator(mat,MATAIJ,&Dmat)); 1891f3b1f45cSBarry Smith 18925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff)); 18935f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN)); 18945f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm)); 18955f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm)); 1896f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1897f3b1f45cSBarry Smith flag = PETSC_FALSE; 1898f3b1f45cSBarry Smith if (v) { 18995f80ce2aSJacob Faibussowitsch CHKERRQ(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))); 19005f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view")); 19015f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view")); 19025f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view")); 1903f3b1f45cSBarry Smith } 1904f3b1f45cSBarry Smith } else if (v) { 19055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n")); 1906f3b1f45cSBarry Smith } 1907f3b1f45cSBarry Smith if (flg) *flg = flag; 19085f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Ddiff)); 19095f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mf)); 19105f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Dmf)); 19115f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Dmat)); 1912f3b1f45cSBarry Smith PetscFunctionReturn(0); 1913f3b1f45cSBarry Smith } 1914f3b1f45cSBarry Smith 1915f3b1f45cSBarry Smith /*@C 19167301b172SPierre Jolivet MatShellTestMultTranspose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function. 1917f3b1f45cSBarry Smith 1918f3b1f45cSBarry Smith Logically Collective on Mat 1919f3b1f45cSBarry Smith 1920f3b1f45cSBarry Smith Input Parameters: 1921f3b1f45cSBarry Smith + mat - the shell matrix 1922f3b1f45cSBarry Smith . f - the function 1923f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1924f3b1f45cSBarry Smith - ctx - an optional context for the function 1925f3b1f45cSBarry Smith 1926f3b1f45cSBarry Smith Output Parameter: 1927f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1928f3b1f45cSBarry Smith 1929f3b1f45cSBarry Smith Options Database: 1930f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1931f3b1f45cSBarry Smith 1932f3b1f45cSBarry Smith Level: advanced 1933f3b1f45cSBarry Smith 193495452b02SPatrick Sanan Fortran Notes: 193595452b02SPatrick Sanan Not supported from Fortran 1936f3b1f45cSBarry Smith 1937f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult() 1938f3b1f45cSBarry Smith @*/ 1939f3b1f45cSBarry Smith PetscErrorCode MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1940f3b1f45cSBarry Smith { 1941f3b1f45cSBarry Smith Vec x,y,z; 1942f3b1f45cSBarry Smith PetscInt m,n,M,N; 1943f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1944f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 194574e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1946f3b1f45cSBarry Smith 1947f3b1f45cSBarry Smith PetscFunctionBegin; 1948f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 19495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v)); 19505f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(mat,&x,&y)); 19515f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(y,&z)); 19525f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(mat,&m,&n)); 19535f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(mat,&M,&N)); 19545f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf)); 19555f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetFunction(mf,f,ctx)); 19565f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetBase(mf,base,NULL)); 19575f80ce2aSJacob Faibussowitsch CHKERRQ(MatComputeOperator(mf,MATAIJ,&Dmf)); 19585f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf)); 19595f80ce2aSJacob Faibussowitsch CHKERRQ(MatComputeOperatorTranspose(mat,MATAIJ,&Dmat)); 1960f3b1f45cSBarry Smith 19615f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff)); 19625f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN)); 19635f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm)); 19645f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm)); 1965f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1966f3b1f45cSBarry Smith flag = PETSC_FALSE; 1967f3b1f45cSBarry Smith if (v) { 19685f80ce2aSJacob Faibussowitsch CHKERRQ(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))); 19695f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view")); 19705f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view")); 19715f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view")); 1972f3b1f45cSBarry Smith } 1973f3b1f45cSBarry Smith } else if (v) { 19745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n")); 1975f3b1f45cSBarry Smith } 1976f3b1f45cSBarry Smith if (flg) *flg = flag; 19775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mf)); 19785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Dmat)); 19795f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Ddiff)); 19805f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Dmf)); 19815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 19825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 19835f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&z)); 1984f3b1f45cSBarry Smith PetscFunctionReturn(0); 1985f3b1f45cSBarry Smith } 1986f3b1f45cSBarry Smith 1987f3b1f45cSBarry Smith /*@C 19880c0fd78eSBarry Smith MatShellSetOperation - Allows user to set a matrix operation for a shell matrix. 1989e51e0e81SBarry Smith 19903f9fe445SBarry Smith Logically Collective on Mat 1991fee21e36SBarry Smith 1992c7fcc2eaSBarry Smith Input Parameters: 1993c7fcc2eaSBarry Smith + mat - the shell matrix 1994c7fcc2eaSBarry Smith . op - the name of the operation 1995789d8953SBarry Smith - g - the function that provides the operation. 1996c7fcc2eaSBarry Smith 199715091d37SBarry Smith Level: advanced 199815091d37SBarry Smith 1999fae171e0SBarry Smith Usage: 2000e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 2001b77ba244SStefano Zampini $ MatCreateShell(comm,m,n,M,N,ctx,&A); 2002b77ba244SStefano Zampini $ MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 20030b627109SLois Curfman McInnes 2004a62d957aSLois Curfman McInnes Notes: 2005e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 20061c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 2007a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 20081c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 2009a62d957aSLois Curfman McInnes 2010a4d85fd7SBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 2011deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 2012deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 2013deebb3c3SLois Curfman McInnes routines, e.g., 2014a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 2015a62d957aSLois Curfman McInnes 20164aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 20174aa34b0aSBarry Smith nonzero on failure. 20184aa34b0aSBarry Smith 2019a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 2020a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 2021a62d957aSLois Curfman McInnes set by MatCreateShell(). 2022a62d957aSLois Curfman McInnes 2023b77ba244SStefano Zampini Use MatSetOperation() to set an operation for any matrix type. For matrix product operations (i.e. MatMat, MatTransposeMat etc) use MatShellSetMatProductOperation() 2024b77ba244SStefano Zampini 202595452b02SPatrick Sanan Fortran Notes: 202695452b02SPatrick Sanan For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 2027c4762a1bSJed Brown generate a matrix. See src/mat/tests/ex120f.F 2028501d9185SBarry Smith 2029b77ba244SStefano Zampini .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts(), MatShellSetMatProductOperation() 2030e51e0e81SBarry Smith @*/ 2031789d8953SBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*g)(void)) 2032e51e0e81SBarry Smith { 20333a40ed3dSBarry Smith PetscFunctionBegin; 20340700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 20355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(mat,"MatShellSetOperation_C",(Mat,MatOperation,void (*)(void)),(mat,op,g))); 20363a40ed3dSBarry Smith PetscFunctionReturn(0); 2037e51e0e81SBarry Smith } 2038f0479e8cSBarry Smith 2039d4bb536fSBarry Smith /*@C 2040d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 2041d4bb536fSBarry Smith 2042c7fcc2eaSBarry Smith Not Collective 2043c7fcc2eaSBarry Smith 2044d4bb536fSBarry Smith Input Parameters: 2045c7fcc2eaSBarry Smith + mat - the shell matrix 2046c7fcc2eaSBarry Smith - op - the name of the operation 2047d4bb536fSBarry Smith 2048d4bb536fSBarry Smith Output Parameter: 2049789d8953SBarry Smith . g - the function that provides the operation. 2050d4bb536fSBarry Smith 205115091d37SBarry Smith Level: advanced 205215091d37SBarry Smith 2053d4bb536fSBarry Smith Notes: 2054e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 2055d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 2056d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 2057d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 2058d4bb536fSBarry Smith 2059d4bb536fSBarry Smith All user-provided functions have the same calling 2060d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 2061d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 2062d4bb536fSBarry Smith routines, e.g., 2063d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 2064d4bb536fSBarry Smith 2065d4bb536fSBarry Smith Within each user-defined routine, the user should call 2066d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 2067d4bb536fSBarry Smith set by MatCreateShell(). 2068d4bb536fSBarry Smith 2069ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 2070d4bb536fSBarry Smith @*/ 2071789d8953SBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**g)(void)) 2072d4bb536fSBarry Smith { 20733a40ed3dSBarry Smith PetscFunctionBegin; 20740700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 20755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscUseMethod(mat,"MatShellGetOperation_C",(Mat,MatOperation,void (**)(void)),(mat,op,g))); 20763a40ed3dSBarry Smith PetscFunctionReturn(0); 2077d4bb536fSBarry Smith } 2078b77ba244SStefano Zampini 2079b77ba244SStefano Zampini /*@ 2080b77ba244SStefano Zampini MatIsShell - Inquires if a matrix is derived from MATSHELL 2081b77ba244SStefano Zampini 2082b77ba244SStefano Zampini Input Parameter: 2083b77ba244SStefano Zampini . mat - the matrix 2084b77ba244SStefano Zampini 2085b77ba244SStefano Zampini Output Parameter: 2086b77ba244SStefano Zampini . flg - the boolean value 2087b77ba244SStefano Zampini 2088b77ba244SStefano Zampini Level: developer 2089b77ba244SStefano Zampini 2090b77ba244SStefano Zampini Notes: in the future, we should allow the object type name to be changed still using the MatShell data structure for other matrices (i.e. MATTRANSPOSEMAT, MATSCHURCOMPLEMENT etc) 2091b77ba244SStefano Zampini 2092b77ba244SStefano Zampini .seealso: MatCreateShell() 2093b77ba244SStefano Zampini @*/ 2094b77ba244SStefano Zampini PetscErrorCode MatIsShell(Mat mat, PetscBool *flg) 2095b77ba244SStefano Zampini { 2096b77ba244SStefano Zampini PetscFunctionBegin; 2097b77ba244SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2098b77ba244SStefano Zampini PetscValidPointer(flg,2); 2099b77ba244SStefano Zampini *flg = (PetscBool)(mat->ops->destroy == MatDestroy_Shell); 2100b77ba244SStefano Zampini PetscFunctionReturn(0); 2101b77ba244SStefano Zampini } 2102