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 */ 62*800f99ffSJeremy L Thompson PetscContainer ctxcontainer; 6388cf3e7dSBarry Smith } Mat_Shell; 640c0fd78eSBarry Smith 6545960306SStefano Zampini /* 6645960306SStefano Zampini Store and scale values on zeroed rows 6745960306SStefano Zampini xx = [x_1, 0], 0 on zeroed columns 6845960306SStefano Zampini */ 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) { 769566063dSJacob Faibussowitsch PetscCall(VecSet(shell->zvals_w,0.0)); 779566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 789566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 799566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals)); 8045960306SStefano Zampini } 8145960306SStefano Zampini if (shell->zcols) { 8245960306SStefano Zampini if (!shell->right_work) { 839566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&shell->right_work,NULL)); 8445960306SStefano Zampini } 859566063dSJacob Faibussowitsch PetscCall(VecCopy(x,shell->right_work)); 869566063dSJacob Faibussowitsch PetscCall(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) { 999566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE)); 1009566063dSJacob Faibussowitsch PetscCall(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) { 1199566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,NULL,&shell->left_work)); 12045960306SStefano Zampini } 1219566063dSJacob Faibussowitsch PetscCall(VecCopy(x,shell->left_work)); 1229566063dSJacob Faibussowitsch PetscCall(VecSet(shell->zvals_w,0.0)); 1239566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE)); 1249566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE)); 1259566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 1269566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 1279566063dSJacob Faibussowitsch PetscCall(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; 1391baa6e33SBarry Smith if (shell->zcols) PetscCall(VecISSet(x,shell->zcols,0.0)); 14045960306SStefano Zampini if (shell->zrows) { 1419566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE)); 1429566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE)); 14345960306SStefano Zampini } 14445960306SStefano Zampini PetscFunctionReturn(0); 14545960306SStefano Zampini } 14645960306SStefano Zampini 1478fe8eb27SJed Brown /* 1480c0fd78eSBarry Smith xx = diag(left)*x 1498fe8eb27SJed Brown */ 1508fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 1518fe8eb27SJed Brown { 1528fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1538fe8eb27SJed Brown 1548fe8eb27SJed Brown PetscFunctionBegin; 1550298fd71SBarry Smith *xx = NULL; 1568fe8eb27SJed Brown if (!shell->left) { 1578fe8eb27SJed Brown *xx = x; 1588fe8eb27SJed Brown } else { 1599566063dSJacob Faibussowitsch if (!shell->left_work) PetscCall(VecDuplicate(shell->left,&shell->left_work)); 1609566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->left_work,x,shell->left)); 1618fe8eb27SJed Brown *xx = shell->left_work; 1628fe8eb27SJed Brown } 1638fe8eb27SJed Brown PetscFunctionReturn(0); 1648fe8eb27SJed Brown } 1658fe8eb27SJed Brown 1660c0fd78eSBarry Smith /* 1670c0fd78eSBarry Smith xx = diag(right)*x 1680c0fd78eSBarry Smith */ 1698fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 1708fe8eb27SJed Brown { 1718fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1728fe8eb27SJed Brown 1738fe8eb27SJed Brown PetscFunctionBegin; 1740298fd71SBarry Smith *xx = NULL; 1758fe8eb27SJed Brown if (!shell->right) { 1768fe8eb27SJed Brown *xx = x; 1778fe8eb27SJed Brown } else { 1789566063dSJacob Faibussowitsch if (!shell->right_work) PetscCall(VecDuplicate(shell->right,&shell->right_work)); 1799566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->right_work,x,shell->right)); 1808fe8eb27SJed Brown *xx = shell->right_work; 1818fe8eb27SJed Brown } 1828fe8eb27SJed Brown PetscFunctionReturn(0); 1838fe8eb27SJed Brown } 1848fe8eb27SJed Brown 1850c0fd78eSBarry Smith /* 1860c0fd78eSBarry Smith x = diag(left)*x 1870c0fd78eSBarry Smith */ 1888fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 1898fe8eb27SJed Brown { 1908fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1918fe8eb27SJed Brown 1928fe8eb27SJed Brown PetscFunctionBegin; 1939566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseMult(x,x,shell->left)); 1948fe8eb27SJed Brown PetscFunctionReturn(0); 1958fe8eb27SJed Brown } 1968fe8eb27SJed Brown 1970c0fd78eSBarry Smith /* 1980c0fd78eSBarry Smith x = diag(right)*x 1990c0fd78eSBarry Smith */ 2008fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 2018fe8eb27SJed Brown { 2028fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2038fe8eb27SJed Brown 2048fe8eb27SJed Brown PetscFunctionBegin; 2059566063dSJacob Faibussowitsch if (shell->right) PetscCall(VecPointwiseMult(x,x,shell->right)); 2068fe8eb27SJed Brown PetscFunctionReturn(0); 2078fe8eb27SJed Brown } 2088fe8eb27SJed Brown 2090c0fd78eSBarry Smith /* 2100c0fd78eSBarry Smith Y = vscale*Y + diag(dshift)*X + vshift*X 2110c0fd78eSBarry Smith 2120c0fd78eSBarry Smith On input Y already contains A*x 2130c0fd78eSBarry Smith */ 2148fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 2158fe8eb27SJed Brown { 2168fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2178fe8eb27SJed Brown 2188fe8eb27SJed Brown PetscFunctionBegin; 2198fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 2208fe8eb27SJed Brown PetscInt i,m; 2218fe8eb27SJed Brown const PetscScalar *x,*d; 2228fe8eb27SJed Brown PetscScalar *y; 2239566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X,&m)); 2249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(shell->dshift,&d)); 2259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 2269566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y,&y)); 2278fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 2289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(shell->dshift,&d)); 2299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 2309566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y,&y)); 2310c0fd78eSBarry Smith } else { 2329566063dSJacob Faibussowitsch PetscCall(VecScale(Y,shell->vscale)); 2338fe8eb27SJed Brown } 2349566063dSJacob Faibussowitsch if (shell->vshift != 0.0) PetscCall(VecAXPY(Y,shell->vshift,X)); /* if test is for non-square matrices */ 2358fe8eb27SJed Brown PetscFunctionReturn(0); 2368fe8eb27SJed Brown } 237e51e0e81SBarry Smith 238*800f99ffSJeremy L Thompson static PetscErrorCode MatShellGetContext_Shell(Mat mat,void *ctx) 239789d8953SBarry Smith { 240*800f99ffSJeremy L Thompson Mat_Shell *shell = (Mat_Shell*)mat->data; 241*800f99ffSJeremy L Thompson 242789d8953SBarry Smith PetscFunctionBegin; 243*800f99ffSJeremy L Thompson if (shell->ctxcontainer) PetscCall(PetscContainerGetPointer(shell->ctxcontainer,(void**)ctx)); 244*800f99ffSJeremy L Thompson else *(void**)ctx = NULL; 245789d8953SBarry Smith PetscFunctionReturn(0); 246789d8953SBarry Smith } 247789d8953SBarry Smith 2489d225801SJed Brown /*@ 249a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 250b4fd4287SBarry Smith 251c7fcc2eaSBarry Smith Not Collective 252c7fcc2eaSBarry Smith 253b4fd4287SBarry Smith Input Parameter: 254b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 255b4fd4287SBarry Smith 256b4fd4287SBarry Smith Output Parameter: 257b4fd4287SBarry Smith . ctx - the user provided context 258b4fd4287SBarry Smith 25915091d37SBarry Smith Level: advanced 26015091d37SBarry Smith 26195452b02SPatrick Sanan Fortran Notes: 26295452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 263daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 264a62d957aSLois Curfman McInnes 265db781477SPatrick Sanan .seealso: `MatCreateShell()`, `MatShellSetOperation()`, `MatShellSetContext()` 2660bc0a719SBarry Smith @*/ 2678fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 268b4fd4287SBarry Smith { 2693a40ed3dSBarry Smith PetscFunctionBegin; 2700700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2714482741eSBarry Smith PetscValidPointer(ctx,2); 272cac4c232SBarry Smith PetscUseMethod(mat,"MatShellGetContext_C",(Mat,void*),(mat,ctx)); 2733a40ed3dSBarry Smith PetscFunctionReturn(0); 274b4fd4287SBarry Smith } 275b4fd4287SBarry Smith 27645960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Local_Shell(Mat mat,PetscInt nr,PetscInt rows[],PetscInt nc,PetscInt cols[],PetscScalar diag,PetscBool rc) 27745960306SStefano Zampini { 27845960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)mat->data; 27945960306SStefano Zampini Vec x = NULL,b = NULL; 28045960306SStefano Zampini IS is1, is2; 28145960306SStefano Zampini const PetscInt *ridxs; 28245960306SStefano Zampini PetscInt *idxs,*gidxs; 28345960306SStefano Zampini PetscInt cum,rst,cst,i; 28445960306SStefano Zampini 28545960306SStefano Zampini PetscFunctionBegin; 28645960306SStefano Zampini if (!shell->zvals) { 2879566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat,NULL,&shell->zvals)); 28845960306SStefano Zampini } 28945960306SStefano Zampini if (!shell->zvals_w) { 2909566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->zvals,&shell->zvals_w)); 29145960306SStefano Zampini } 2929566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat,&rst,NULL)); 2939566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(mat,&cst,NULL)); 29445960306SStefano Zampini 29545960306SStefano Zampini /* Expand/create index set of zeroed rows */ 2969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nr,&idxs)); 29745960306SStefano Zampini for (i = 0; i < nr; i++) idxs[i] = rows[i] + rst; 2989566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nr,idxs,PETSC_OWN_POINTER,&is1)); 2999566063dSJacob Faibussowitsch PetscCall(ISSort(is1)); 3009566063dSJacob Faibussowitsch PetscCall(VecISSet(shell->zvals,is1,diag)); 30145960306SStefano Zampini if (shell->zrows) { 3029566063dSJacob Faibussowitsch PetscCall(ISSum(shell->zrows,is1,&is2)); 3039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zrows)); 3049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 30545960306SStefano Zampini shell->zrows = is2; 30645960306SStefano Zampini } else shell->zrows = is1; 30745960306SStefano Zampini 30845960306SStefano Zampini /* Create scatters for diagonal values communications */ 3099566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_c)); 3109566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_r)); 31145960306SStefano Zampini 31245960306SStefano Zampini /* row scatter: from/to left vector */ 3139566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat,&x,&b)); 3149566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(b,shell->zrows,shell->zvals_w,shell->zrows,&shell->zvals_sct_r)); 31545960306SStefano Zampini 31645960306SStefano Zampini /* col scatter: from right vector to left vector */ 3179566063dSJacob Faibussowitsch PetscCall(ISGetIndices(shell->zrows,&ridxs)); 3189566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(shell->zrows,&nr)); 3199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nr,&gidxs)); 32045960306SStefano Zampini for (i = 0, cum = 0; i < nr; i++) { 32145960306SStefano Zampini if (ridxs[i] >= mat->cmap->N) continue; 32245960306SStefano Zampini gidxs[cum] = ridxs[i]; 32345960306SStefano Zampini cum++; 32445960306SStefano Zampini } 3259566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)mat),cum,gidxs,PETSC_OWN_POINTER,&is1)); 3269566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,is1,shell->zvals_w,is1,&shell->zvals_sct_c)); 3279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 3289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 33045960306SStefano Zampini 33145960306SStefano Zampini /* Expand/create index set of zeroed columns */ 33245960306SStefano Zampini if (rc) { 3339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nc,&idxs)); 33445960306SStefano Zampini for (i = 0; i < nc; i++) idxs[i] = cols[i] + cst; 3359566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nc,idxs,PETSC_OWN_POINTER,&is1)); 3369566063dSJacob Faibussowitsch PetscCall(ISSort(is1)); 33745960306SStefano Zampini if (shell->zcols) { 3389566063dSJacob Faibussowitsch PetscCall(ISSum(shell->zcols,is1,&is2)); 3399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zcols)); 3409566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 34145960306SStefano Zampini shell->zcols = is2; 34245960306SStefano Zampini } else shell->zcols = is1; 34345960306SStefano Zampini } 34445960306SStefano Zampini PetscFunctionReturn(0); 34545960306SStefano Zampini } 34645960306SStefano Zampini 34745960306SStefano Zampini static PetscErrorCode MatZeroRows_Shell(Mat mat,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 34845960306SStefano Zampini { 349ef5c7bd2SStefano Zampini Mat_Shell *shell = (Mat_Shell*)mat->data; 35045960306SStefano Zampini PetscInt nr, *lrows; 35145960306SStefano Zampini 35245960306SStefano Zampini PetscFunctionBegin; 35345960306SStefano Zampini if (x && b) { 35445960306SStefano Zampini Vec xt; 35545960306SStefano Zampini PetscScalar *vals; 35645960306SStefano Zampini PetscInt *gcols,i,st,nl,nc; 35745960306SStefano Zampini 3589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&gcols)); 35945960306SStefano Zampini for (i = 0, nc = 0; i < n; i++) if (rows[i] < mat->cmap->N) gcols[nc++] = rows[i]; 36045960306SStefano Zampini 3619566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat,&xt,NULL)); 3629566063dSJacob Faibussowitsch PetscCall(VecCopy(x,xt)); 3639566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nc,&vals)); 3649566063dSJacob Faibussowitsch PetscCall(VecSetValues(xt,nc,gcols,vals,INSERT_VALUES)); /* xt = [x1, 0] */ 3659566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 3669566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(xt)); 3679566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(xt)); 3689566063dSJacob Faibussowitsch PetscCall(VecAYPX(xt,-1.0,x)); /* xt = [0, x2] */ 36945960306SStefano Zampini 3709566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(xt,&st,NULL)); 3719566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xt,&nl)); 3729566063dSJacob Faibussowitsch PetscCall(VecGetArray(xt,&vals)); 37345960306SStefano Zampini for (i = 0; i < nl; i++) { 37445960306SStefano Zampini PetscInt g = i + st; 37545960306SStefano Zampini if (g > mat->rmap->N) continue; 37645960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 3779566063dSJacob Faibussowitsch PetscCall(VecSetValue(b,g,diag*vals[i],INSERT_VALUES)); 37845960306SStefano Zampini } 3799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xt,&vals)); 3809566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b)); 3819566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b)); /* b = [b1, x2 * diag] */ 3829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xt)); 3839566063dSJacob Faibussowitsch PetscCall(PetscFree(gcols)); 38445960306SStefano Zampini } 3859566063dSJacob Faibussowitsch PetscCall(PetscLayoutMapLocal(mat->rmap,n,rows,&nr,&lrows,NULL)); 3869566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns_Local_Shell(mat,nr,lrows,0,NULL,diag,PETSC_FALSE)); 3871baa6e33SBarry Smith if (shell->axpy) PetscCall(MatZeroRows(shell->axpy,n,rows,0.0,NULL,NULL)); 3889566063dSJacob Faibussowitsch PetscCall(PetscFree(lrows)); 38945960306SStefano Zampini PetscFunctionReturn(0); 39045960306SStefano Zampini } 39145960306SStefano Zampini 39245960306SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Shell(Mat mat,PetscInt n,const PetscInt rowscols[],PetscScalar diag,Vec x,Vec b) 39345960306SStefano Zampini { 394ef5c7bd2SStefano Zampini Mat_Shell *shell = (Mat_Shell*)mat->data; 39545960306SStefano Zampini PetscInt *lrows, *lcols; 39645960306SStefano Zampini PetscInt nr, nc; 39745960306SStefano Zampini PetscBool congruent; 39845960306SStefano Zampini 39945960306SStefano Zampini PetscFunctionBegin; 40045960306SStefano Zampini if (x && b) { 40145960306SStefano Zampini Vec xt, bt; 40245960306SStefano Zampini PetscScalar *vals; 40345960306SStefano Zampini PetscInt *grows,*gcols,i,st,nl; 40445960306SStefano Zampini 4059566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n,&grows,n,&gcols)); 40645960306SStefano Zampini for (i = 0, nr = 0; i < n; i++) if (rowscols[i] < mat->rmap->N) grows[nr++] = rowscols[i]; 40745960306SStefano Zampini for (i = 0, nc = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) gcols[nc++] = rowscols[i]; 4089566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(n,&vals)); 40945960306SStefano Zampini 4109566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat,&xt,&bt)); 4119566063dSJacob Faibussowitsch PetscCall(VecCopy(x,xt)); 4129566063dSJacob Faibussowitsch PetscCall(VecSetValues(xt,nc,gcols,vals,INSERT_VALUES)); /* xt = [x1, 0] */ 4139566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(xt)); 4149566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(xt)); 4159566063dSJacob Faibussowitsch PetscCall(VecAXPY(xt,-1.0,x)); /* xt = [0, -x2] */ 4169566063dSJacob Faibussowitsch PetscCall(MatMult(mat,xt,bt)); /* bt = [-A12*x2,-A22*x2] */ 4179566063dSJacob Faibussowitsch PetscCall(VecSetValues(bt,nr,grows,vals,INSERT_VALUES)); /* bt = [-A12*x2,0] */ 4189566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bt)); 4199566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bt)); 4209566063dSJacob Faibussowitsch PetscCall(VecAXPY(b,1.0,bt)); /* b = [b1 - A12*x2, b2] */ 4219566063dSJacob Faibussowitsch PetscCall(VecSetValues(bt,nr,grows,vals,INSERT_VALUES)); /* b = [b1 - A12*x2, 0] */ 4229566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bt)); 4239566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bt)); 4249566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 42545960306SStefano Zampini 4269566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(xt,&st,NULL)); 4279566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xt,&nl)); 4289566063dSJacob Faibussowitsch PetscCall(VecGetArray(xt,&vals)); 42945960306SStefano Zampini for (i = 0; i < nl; i++) { 43045960306SStefano Zampini PetscInt g = i + st; 43145960306SStefano Zampini if (g > mat->rmap->N) continue; 43245960306SStefano Zampini if (PetscAbsScalar(vals[i]) == 0.0) continue; 4339566063dSJacob Faibussowitsch PetscCall(VecSetValue(b,g,-diag*vals[i],INSERT_VALUES)); 43445960306SStefano Zampini } 4359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xt,&vals)); 4369566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b)); 4379566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b)); /* b = [b1 - A12*x2, x2 * diag] */ 4389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xt)); 4399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bt)); 4409566063dSJacob Faibussowitsch PetscCall(PetscFree2(grows,gcols)); 44145960306SStefano Zampini } 4429566063dSJacob Faibussowitsch PetscCall(PetscLayoutMapLocal(mat->rmap,n,rowscols,&nr,&lrows,NULL)); 4439566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(mat,&congruent)); 44445960306SStefano Zampini if (congruent) { 44545960306SStefano Zampini nc = nr; 44645960306SStefano Zampini lcols = lrows; 44745960306SStefano Zampini } else { /* MatZeroRowsColumns implicitly assumes the rowscols indices are for a square matrix, here we handle a more general case */ 44845960306SStefano Zampini PetscInt i,nt,*t; 44945960306SStefano Zampini 4509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&t)); 45145960306SStefano Zampini for (i = 0, nt = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) t[nt++] = rowscols[i]; 4529566063dSJacob Faibussowitsch PetscCall(PetscLayoutMapLocal(mat->cmap,nt,t,&nc,&lcols,NULL)); 4539566063dSJacob Faibussowitsch PetscCall(PetscFree(t)); 45445960306SStefano Zampini } 4559566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns_Local_Shell(mat,nr,lrows,nc,lcols,diag,PETSC_TRUE)); 45645960306SStefano Zampini if (!congruent) { 4579566063dSJacob Faibussowitsch PetscCall(PetscFree(lcols)); 45845960306SStefano Zampini } 4599566063dSJacob Faibussowitsch PetscCall(PetscFree(lrows)); 4601baa6e33SBarry Smith if (shell->axpy) PetscCall(MatZeroRowsColumns(shell->axpy,n,rowscols,0.0,NULL,NULL)); 46145960306SStefano Zampini PetscFunctionReturn(0); 46245960306SStefano Zampini } 46345960306SStefano Zampini 464*800f99ffSJeremy L Thompson static PetscErrorCode MatDestroy_Shell(Mat mat) 465e51e0e81SBarry Smith { 466bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 467b77ba244SStefano Zampini MatShellMatFunctionList matmat; 468ed3cc1f0SBarry Smith 4693a40ed3dSBarry Smith PetscFunctionBegin; 4701baa6e33SBarry Smith if (shell->ops->destroy) PetscCall((*shell->ops->destroy)(mat)); 4719566063dSJacob Faibussowitsch PetscCall(PetscMemzero(shell->ops,sizeof(struct _MatShellOps))); 4729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->left)); 4739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->right)); 4749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->dshift)); 4759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->left_work)); 4769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->right_work)); 4779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->left_add_work)); 4789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->right_add_work)); 4799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->axpy_left)); 4809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->axpy_right)); 4819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&shell->axpy)); 4829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->zvals_w)); 4839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->zvals)); 4849566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_c)); 4859566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_r)); 4869566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zrows)); 4879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zcols)); 488b77ba244SStefano Zampini 489b77ba244SStefano Zampini matmat = shell->matmat; 490b77ba244SStefano Zampini while (matmat) { 491b77ba244SStefano Zampini MatShellMatFunctionList next = matmat->next; 492b77ba244SStefano Zampini 4939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,matmat->composedname,NULL)); 4949566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat->composedname)); 4959566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat->resultname)); 4969566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat)); 497b77ba244SStefano Zampini matmat = next; 498b77ba244SStefano Zampini } 499*800f99ffSJeremy L Thompson PetscCall(MatShellSetContext(mat,NULL)); 5009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatShellGetContext_C",NULL)); 5019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetContext_C",NULL)); 502*800f99ffSJeremy L Thompson PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetContextDestroy_C",NULL)); 5039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetVecType_C",NULL)); 5049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetManageScalingShifts_C",NULL)); 5059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetOperation_C",NULL)); 5069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatShellGetOperation_C",NULL)); 5079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatShellSetMatProductOperation_C",NULL)); 5089566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 509b77ba244SStefano Zampini PetscFunctionReturn(0); 510b77ba244SStefano Zampini } 511b77ba244SStefano Zampini 512b77ba244SStefano Zampini typedef struct { 513b77ba244SStefano Zampini PetscErrorCode (*numeric)(Mat,Mat,Mat,void*); 514b77ba244SStefano Zampini PetscErrorCode (*destroy)(void*); 515b77ba244SStefano Zampini void *userdata; 516b77ba244SStefano Zampini Mat B; 517b77ba244SStefano Zampini Mat Bt; 518b77ba244SStefano Zampini Mat axpy; 519b77ba244SStefano Zampini } MatMatDataShell; 520b77ba244SStefano Zampini 521b77ba244SStefano Zampini static PetscErrorCode DestroyMatMatDataShell(void *data) 522b77ba244SStefano Zampini { 523b77ba244SStefano Zampini MatMatDataShell *mmdata = (MatMatDataShell *)data; 524b77ba244SStefano Zampini 525b77ba244SStefano Zampini PetscFunctionBegin; 5261baa6e33SBarry Smith if (mmdata->destroy) PetscCall((*mmdata->destroy)(mmdata->userdata)); 5279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mmdata->B)); 5289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mmdata->Bt)); 5299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mmdata->axpy)); 5309566063dSJacob Faibussowitsch PetscCall(PetscFree(mmdata)); 531b77ba244SStefano Zampini PetscFunctionReturn(0); 532b77ba244SStefano Zampini } 533b77ba244SStefano Zampini 534b77ba244SStefano Zampini static PetscErrorCode MatProductNumeric_Shell_X(Mat D) 535b77ba244SStefano Zampini { 536b77ba244SStefano Zampini Mat_Product *product; 537b77ba244SStefano Zampini Mat A, B; 538b77ba244SStefano Zampini MatMatDataShell *mdata; 539b77ba244SStefano Zampini PetscScalar zero = 0.0; 540b77ba244SStefano Zampini 541b77ba244SStefano Zampini PetscFunctionBegin; 542b77ba244SStefano Zampini MatCheckProduct(D,1); 543b77ba244SStefano Zampini product = D->product; 54428b400f6SJacob Faibussowitsch PetscCheck(product->data,PetscObjectComm((PetscObject)D),PETSC_ERR_PLIB,"Product data empty"); 545b77ba244SStefano Zampini A = product->A; 546b77ba244SStefano Zampini B = product->B; 547b77ba244SStefano Zampini mdata = (MatMatDataShell*)product->data; 548b77ba244SStefano Zampini if (mdata->numeric) { 549b77ba244SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 550b77ba244SStefano Zampini PetscErrorCode (*stashsym)(Mat) = D->ops->productsymbolic; 551b77ba244SStefano Zampini PetscErrorCode (*stashnum)(Mat) = D->ops->productnumeric; 552b77ba244SStefano Zampini PetscBool useBmdata = PETSC_FALSE, newB = PETSC_TRUE; 553b77ba244SStefano Zampini 554b77ba244SStefano Zampini if (shell->managescalingshifts) { 55508401ef6SPierre Jolivet PetscCheck(!shell->zcols && !shell->zrows,PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProduct not supported with zeroed rows/columns"); 556b77ba244SStefano Zampini if (shell->right || shell->left) { 557b77ba244SStefano Zampini useBmdata = PETSC_TRUE; 558b77ba244SStefano Zampini if (!mdata->B) { 5599566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_SHARE_NONZERO_PATTERN,&mdata->B)); 560b77ba244SStefano Zampini } else { 561b77ba244SStefano Zampini newB = PETSC_FALSE; 562b77ba244SStefano Zampini } 5639566063dSJacob Faibussowitsch PetscCall(MatCopy(B,mdata->B,SAME_NONZERO_PATTERN)); 564b77ba244SStefano Zampini } 565b77ba244SStefano Zampini switch (product->type) { 566b77ba244SStefano Zampini case MATPRODUCT_AB: /* s L A R B + v L R B + L D R B */ 5671baa6e33SBarry Smith if (shell->right) PetscCall(MatDiagonalScale(mdata->B,shell->right,NULL)); 568b77ba244SStefano Zampini break; 569b77ba244SStefano Zampini case MATPRODUCT_AtB: /* s R A^t L B + v R L B + R D L B */ 5701baa6e33SBarry Smith if (shell->left) PetscCall(MatDiagonalScale(mdata->B,shell->left,NULL)); 571b77ba244SStefano Zampini break; 572b77ba244SStefano Zampini case MATPRODUCT_ABt: /* s L A R B^t + v L R B^t + L D R B^t */ 5731baa6e33SBarry Smith if (shell->right) PetscCall(MatDiagonalScale(mdata->B,NULL,shell->right)); 574b77ba244SStefano Zampini break; 575b77ba244SStefano Zampini case MATPRODUCT_RARt: /* s B L A R B^t + v B L R B^t + B L D R B^t */ 576b77ba244SStefano Zampini if (shell->right && shell->left) { 577b77ba244SStefano Zampini PetscBool flg; 578b77ba244SStefano Zampini 5799566063dSJacob Faibussowitsch PetscCall(VecEqual(shell->right,shell->left,&flg)); 58028b400f6SJacob 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); 581b77ba244SStefano Zampini } 5821baa6e33SBarry Smith if (shell->right) PetscCall(MatDiagonalScale(mdata->B,NULL,shell->right)); 583b77ba244SStefano Zampini break; 584b77ba244SStefano Zampini case MATPRODUCT_PtAP: /* s B^t L A R B + v B^t L R B + B^t L D R B */ 585b77ba244SStefano Zampini if (shell->right && shell->left) { 586b77ba244SStefano Zampini PetscBool flg; 587b77ba244SStefano Zampini 5889566063dSJacob Faibussowitsch PetscCall(VecEqual(shell->right,shell->left,&flg)); 58928b400f6SJacob 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); 590b77ba244SStefano Zampini } 5911baa6e33SBarry Smith if (shell->right) PetscCall(MatDiagonalScale(mdata->B,shell->right,NULL)); 592b77ba244SStefano Zampini break; 59398921bdaSJacob 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); 594b77ba244SStefano Zampini } 595b77ba244SStefano Zampini } 596b77ba244SStefano Zampini /* allow the user to call MatMat operations on D */ 597b77ba244SStefano Zampini D->product = NULL; 598b77ba244SStefano Zampini D->ops->productsymbolic = NULL; 599b77ba244SStefano Zampini D->ops->productnumeric = NULL; 600b77ba244SStefano Zampini 6019566063dSJacob Faibussowitsch PetscCall((*mdata->numeric)(A,useBmdata ? mdata->B : B,D,mdata->userdata)); 602b77ba244SStefano Zampini 603b77ba244SStefano Zampini /* clear any leftover user data and restore D pointers */ 6049566063dSJacob Faibussowitsch PetscCall(MatProductClear(D)); 605b77ba244SStefano Zampini D->ops->productsymbolic = stashsym; 606b77ba244SStefano Zampini D->ops->productnumeric = stashnum; 607b77ba244SStefano Zampini D->product = product; 608b77ba244SStefano Zampini 609b77ba244SStefano Zampini if (shell->managescalingshifts) { 6109566063dSJacob Faibussowitsch PetscCall(MatScale(D,shell->vscale)); 611b77ba244SStefano Zampini switch (product->type) { 612b77ba244SStefano Zampini case MATPRODUCT_AB: /* s L A R B + v L R B + L D R B */ 613b77ba244SStefano Zampini case MATPRODUCT_ABt: /* s L A R B^t + v L R B^t + L D R B^t */ 614b77ba244SStefano Zampini if (shell->left) { 6159566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(D,shell->left,NULL)); 616b77ba244SStefano Zampini if (shell->dshift || shell->vshift != zero) { 6179566063dSJacob Faibussowitsch if (!shell->left_work) PetscCall(MatCreateVecs(A,NULL,&shell->left_work)); 618b77ba244SStefano Zampini if (shell->dshift) { 6199566063dSJacob Faibussowitsch PetscCall(VecCopy(shell->dshift,shell->left_work)); 6209566063dSJacob Faibussowitsch PetscCall(VecShift(shell->left_work,shell->vshift)); 6219566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->left_work,shell->left_work,shell->left)); 622b77ba244SStefano Zampini } else { 6239566063dSJacob Faibussowitsch PetscCall(VecSet(shell->left_work,shell->vshift)); 624b77ba244SStefano Zampini } 625b77ba244SStefano Zampini if (product->type == MATPRODUCT_ABt) { 626b77ba244SStefano Zampini MatReuse reuse = mdata->Bt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX; 627b77ba244SStefano Zampini MatStructure str = mdata->Bt ? SUBSET_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 628b77ba244SStefano Zampini 6299566063dSJacob Faibussowitsch PetscCall(MatTranspose(mdata->B,reuse,&mdata->Bt)); 6309566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(mdata->Bt,shell->left_work,NULL)); 6319566063dSJacob Faibussowitsch PetscCall(MatAXPY(D,1.0,mdata->Bt,str)); 632b77ba244SStefano Zampini } else { 633b77ba244SStefano Zampini MatStructure str = newB ? DIFFERENT_NONZERO_PATTERN : SUBSET_NONZERO_PATTERN; 634b77ba244SStefano Zampini 6359566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(mdata->B,shell->left_work,NULL)); 6369566063dSJacob Faibussowitsch PetscCall(MatAXPY(D,1.0,mdata->B,str)); 637b77ba244SStefano Zampini } 638b77ba244SStefano Zampini } 639b77ba244SStefano Zampini } 640b77ba244SStefano Zampini break; 641b77ba244SStefano Zampini case MATPRODUCT_AtB: /* s R A^t L B + v R L B + R D L B */ 642b77ba244SStefano Zampini if (shell->right) { 6439566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(D,shell->right,NULL)); 644b77ba244SStefano Zampini if (shell->dshift || shell->vshift != zero) { 645b77ba244SStefano Zampini MatStructure str = newB ? DIFFERENT_NONZERO_PATTERN : SUBSET_NONZERO_PATTERN; 646b77ba244SStefano Zampini 6479566063dSJacob Faibussowitsch if (!shell->right_work) PetscCall(MatCreateVecs(A,&shell->right_work,NULL)); 648b77ba244SStefano Zampini if (shell->dshift) { 6499566063dSJacob Faibussowitsch PetscCall(VecCopy(shell->dshift,shell->right_work)); 6509566063dSJacob Faibussowitsch PetscCall(VecShift(shell->right_work,shell->vshift)); 6519566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->right_work,shell->right_work,shell->right)); 652b77ba244SStefano Zampini } else { 6539566063dSJacob Faibussowitsch PetscCall(VecSet(shell->right_work,shell->vshift)); 654b77ba244SStefano Zampini } 6559566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(mdata->B,shell->right_work,NULL)); 6569566063dSJacob Faibussowitsch PetscCall(MatAXPY(D,1.0,mdata->B,str)); 657b77ba244SStefano Zampini } 658b77ba244SStefano Zampini } 659b77ba244SStefano Zampini break; 660b77ba244SStefano Zampini case MATPRODUCT_PtAP: /* s B^t L A R B + v B^t L R B + B^t L D R B */ 661b77ba244SStefano Zampini case MATPRODUCT_RARt: /* s B L A R B^t + v B L R B^t + B L D R B^t */ 662aed4548fSBarry Smith PetscCheck(!shell->dshift && shell->vshift == zero,PetscObjectComm((PetscObject)D),PETSC_ERR_SUP,"MatProductSymbolic type %s not supported for %s and %s matrices with diagonal shift",MatProductTypes[product->type],((PetscObject)A)->type_name,((PetscObject)B)->type_name); 663b77ba244SStefano Zampini break; 66498921bdaSJacob 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); 665b77ba244SStefano Zampini } 666b77ba244SStefano Zampini if (shell->axpy && shell->axpy_vscale != zero) { 667b77ba244SStefano Zampini Mat X; 668b77ba244SStefano Zampini PetscObjectState axpy_state; 669b77ba244SStefano Zampini MatStructure str = DIFFERENT_NONZERO_PATTERN; /* not sure it is safe to ever use SUBSET_NONZERO_PATTERN */ 670b77ba244SStefano Zampini 6719566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(shell->axpy,&X)); 6729566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X,&axpy_state)); 67308401ef6SPierre Jolivet PetscCheck(shell->axpy_state == axpy_state,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)"); 674b77ba244SStefano Zampini if (!mdata->axpy) { 675b77ba244SStefano Zampini str = DIFFERENT_NONZERO_PATTERN; 6769566063dSJacob Faibussowitsch PetscCall(MatProductCreate(shell->axpy,B,NULL,&mdata->axpy)); 6779566063dSJacob Faibussowitsch PetscCall(MatProductSetType(mdata->axpy,product->type)); 6789566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(mdata->axpy)); 6799566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(mdata->axpy)); 680b77ba244SStefano Zampini } else { /* May be that shell->axpy has changed */ 681b77ba244SStefano Zampini PetscBool flg; 682b77ba244SStefano Zampini 6839566063dSJacob Faibussowitsch PetscCall(MatProductReplaceMats(shell->axpy,B,NULL,mdata->axpy)); 6849566063dSJacob Faibussowitsch PetscCall(MatHasOperation(mdata->axpy,MATOP_PRODUCTSYMBOLIC,&flg)); 685b77ba244SStefano Zampini if (!flg) { 686b77ba244SStefano Zampini str = DIFFERENT_NONZERO_PATTERN; 6879566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(mdata->axpy)); 6889566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(mdata->axpy)); 689b77ba244SStefano Zampini } 690b77ba244SStefano Zampini } 6919566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(mdata->axpy)); 6929566063dSJacob Faibussowitsch PetscCall(MatAXPY(D,shell->axpy_vscale,mdata->axpy,str)); 693b77ba244SStefano Zampini } 694b77ba244SStefano Zampini } 695b77ba244SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)D),PETSC_ERR_PLIB,"Missing numeric operation"); 696b77ba244SStefano Zampini PetscFunctionReturn(0); 697b77ba244SStefano Zampini } 698b77ba244SStefano Zampini 699b77ba244SStefano Zampini static PetscErrorCode MatProductSymbolic_Shell_X(Mat D) 700b77ba244SStefano Zampini { 701b77ba244SStefano Zampini Mat_Product *product; 702b77ba244SStefano Zampini Mat A,B; 703b77ba244SStefano Zampini MatShellMatFunctionList matmat; 704b77ba244SStefano Zampini Mat_Shell *shell; 705b77ba244SStefano Zampini PetscBool flg; 706b77ba244SStefano Zampini char composedname[256]; 707b77ba244SStefano Zampini MatMatDataShell *mdata; 708b77ba244SStefano Zampini 709b77ba244SStefano Zampini PetscFunctionBegin; 710b77ba244SStefano Zampini MatCheckProduct(D,1); 711b77ba244SStefano Zampini product = D->product; 71228b400f6SJacob Faibussowitsch PetscCheck(!product->data,PetscObjectComm((PetscObject)D),PETSC_ERR_PLIB,"Product data not empty"); 713b77ba244SStefano Zampini A = product->A; 714b77ba244SStefano Zampini B = product->B; 715b77ba244SStefano Zampini shell = (Mat_Shell*)A->data; 716b77ba244SStefano Zampini matmat = shell->matmat; 7179566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(composedname,sizeof(composedname),"MatProductSetFromOptions_%s_%s_C",((PetscObject)A)->type_name,((PetscObject)B)->type_name)); 718b77ba244SStefano Zampini while (matmat) { 7199566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(composedname,matmat->composedname,&flg)); 720b77ba244SStefano Zampini flg = (PetscBool)(flg && (matmat->ptype == product->type)); 721b77ba244SStefano Zampini if (flg) break; 722b77ba244SStefano Zampini matmat = matmat->next; 723b77ba244SStefano Zampini } 72428b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Composedname \"%s\" for product type %s not found",composedname,MatProductTypes[product->type]); 725b77ba244SStefano Zampini switch (product->type) { 726b77ba244SStefano Zampini case MATPRODUCT_AB: 7279566063dSJacob Faibussowitsch PetscCall(MatSetSizes(D,A->rmap->n,B->cmap->n,A->rmap->N,B->cmap->N)); 728b77ba244SStefano Zampini break; 729b77ba244SStefano Zampini case MATPRODUCT_AtB: 7309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(D,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N)); 731b77ba244SStefano Zampini break; 732b77ba244SStefano Zampini case MATPRODUCT_ABt: 7339566063dSJacob Faibussowitsch PetscCall(MatSetSizes(D,A->rmap->n,B->rmap->n,A->rmap->N,B->rmap->N)); 734b77ba244SStefano Zampini break; 735b77ba244SStefano Zampini case MATPRODUCT_RARt: 7369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(D,B->rmap->n,B->rmap->n,B->rmap->N,B->rmap->N)); 737b77ba244SStefano Zampini break; 738b77ba244SStefano Zampini case MATPRODUCT_PtAP: 7399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(D,B->cmap->n,B->cmap->n,B->cmap->N,B->cmap->N)); 740b77ba244SStefano Zampini break; 74198921bdaSJacob 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); 742b77ba244SStefano Zampini } 743b77ba244SStefano Zampini /* respect users who passed in a matrix for which resultname is the base type */ 744b77ba244SStefano Zampini if (matmat->resultname) { 7459566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)D,matmat->resultname,&flg)); 746b77ba244SStefano Zampini if (!flg) { 7479566063dSJacob Faibussowitsch PetscCall(MatSetType(D,matmat->resultname)); 748b77ba244SStefano Zampini } 749b77ba244SStefano Zampini } 750b77ba244SStefano Zampini /* If matrix type was not set or different, we need to reset this pointers */ 751b77ba244SStefano Zampini D->ops->productsymbolic = MatProductSymbolic_Shell_X; 752b77ba244SStefano Zampini D->ops->productnumeric = MatProductNumeric_Shell_X; 753b77ba244SStefano Zampini /* attach product data */ 7549566063dSJacob Faibussowitsch PetscCall(PetscNew(&mdata)); 755b77ba244SStefano Zampini mdata->numeric = matmat->numeric; 756b77ba244SStefano Zampini mdata->destroy = matmat->destroy; 757b77ba244SStefano Zampini if (matmat->symbolic) { 7589566063dSJacob Faibussowitsch PetscCall((*matmat->symbolic)(A,B,D,&mdata->userdata)); 759b77ba244SStefano Zampini } else { /* call general setup if symbolic operation not provided */ 7609566063dSJacob Faibussowitsch PetscCall(MatSetUp(D)); 761b77ba244SStefano Zampini } 76228b400f6SJacob Faibussowitsch PetscCheck(D->product,PetscObjectComm((PetscObject)D),PETSC_ERR_COR,"Product disappeared after user symbolic phase"); 76328b400f6SJacob Faibussowitsch PetscCheck(!D->product->data,PetscObjectComm((PetscObject)D),PETSC_ERR_COR,"Product data not empty after user symbolic phase"); 764b77ba244SStefano Zampini D->product->data = mdata; 765b77ba244SStefano Zampini D->product->destroy = DestroyMatMatDataShell; 766b77ba244SStefano Zampini /* Be sure to reset these pointers if the user did something unexpected */ 767b77ba244SStefano Zampini D->ops->productsymbolic = MatProductSymbolic_Shell_X; 768b77ba244SStefano Zampini D->ops->productnumeric = MatProductNumeric_Shell_X; 769b77ba244SStefano Zampini PetscFunctionReturn(0); 770b77ba244SStefano Zampini } 771b77ba244SStefano Zampini 772b77ba244SStefano Zampini static PetscErrorCode MatProductSetFromOptions_Shell_X(Mat D) 773b77ba244SStefano Zampini { 774b77ba244SStefano Zampini Mat_Product *product; 775b77ba244SStefano Zampini Mat A,B; 776b77ba244SStefano Zampini MatShellMatFunctionList matmat; 777b77ba244SStefano Zampini Mat_Shell *shell; 778b77ba244SStefano Zampini PetscBool flg; 779b77ba244SStefano Zampini char composedname[256]; 780b77ba244SStefano Zampini 781b77ba244SStefano Zampini PetscFunctionBegin; 782b77ba244SStefano Zampini MatCheckProduct(D,1); 783b77ba244SStefano Zampini product = D->product; 784b77ba244SStefano Zampini A = product->A; 785b77ba244SStefano Zampini B = product->B; 7869566063dSJacob Faibussowitsch PetscCall(MatIsShell(A,&flg)); 787b77ba244SStefano Zampini if (!flg) PetscFunctionReturn(0); 788b77ba244SStefano Zampini shell = (Mat_Shell*)A->data; 789b77ba244SStefano Zampini matmat = shell->matmat; 7909566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(composedname,sizeof(composedname),"MatProductSetFromOptions_%s_%s_C",((PetscObject)A)->type_name,((PetscObject)B)->type_name)); 791b77ba244SStefano Zampini while (matmat) { 7929566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(composedname,matmat->composedname,&flg)); 793b77ba244SStefano Zampini flg = (PetscBool)(flg && (matmat->ptype == product->type)); 794b77ba244SStefano Zampini if (flg) break; 795b77ba244SStefano Zampini matmat = matmat->next; 796b77ba244SStefano Zampini } 797b77ba244SStefano Zampini if (flg) { D->ops->productsymbolic = MatProductSymbolic_Shell_X; } 7989566063dSJacob Faibussowitsch else PetscCall(PetscInfo(D," symbolic product %s not registered for product type %s\n",composedname,MatProductTypes[product->type])); 799b77ba244SStefano Zampini PetscFunctionReturn(0); 800b77ba244SStefano Zampini } 801b77ba244SStefano Zampini 802b77ba244SStefano 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) 803b77ba244SStefano Zampini { 804b77ba244SStefano Zampini PetscBool flg; 805b77ba244SStefano Zampini Mat_Shell *shell; 806b77ba244SStefano Zampini MatShellMatFunctionList matmat; 807b77ba244SStefano Zampini 808b77ba244SStefano Zampini PetscFunctionBegin; 80928b400f6SJacob Faibussowitsch PetscCheck(numeric,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing numeric routine"); 81028b400f6SJacob Faibussowitsch PetscCheck(composedname,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing composed name"); 811b77ba244SStefano Zampini 812b77ba244SStefano Zampini /* add product callback */ 813b77ba244SStefano Zampini shell = (Mat_Shell*)A->data; 814b77ba244SStefano Zampini matmat = shell->matmat; 815b77ba244SStefano Zampini if (!matmat) { 8169566063dSJacob Faibussowitsch PetscCall(PetscNew(&shell->matmat)); 817b77ba244SStefano Zampini matmat = shell->matmat; 818b77ba244SStefano Zampini } else { 819b77ba244SStefano Zampini MatShellMatFunctionList entry = matmat; 820b77ba244SStefano Zampini while (entry) { 8219566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(composedname,entry->composedname,&flg)); 822b77ba244SStefano Zampini flg = (PetscBool)(flg && (entry->ptype == ptype)); 823b77ba244SStefano Zampini matmat = entry; 8242e956fe4SStefano Zampini if (flg) goto set; 825b77ba244SStefano Zampini entry = entry->next; 826b77ba244SStefano Zampini } 8279566063dSJacob Faibussowitsch PetscCall(PetscNew(&matmat->next)); 828b77ba244SStefano Zampini matmat = matmat->next; 829b77ba244SStefano Zampini } 830b77ba244SStefano Zampini 831843d480fSStefano Zampini set: 832b77ba244SStefano Zampini matmat->symbolic = symbolic; 833b77ba244SStefano Zampini matmat->numeric = numeric; 834b77ba244SStefano Zampini matmat->destroy = destroy; 835b77ba244SStefano Zampini matmat->ptype = ptype; 8369566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat->composedname)); 8379566063dSJacob Faibussowitsch PetscCall(PetscFree(matmat->resultname)); 8389566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(composedname,&matmat->composedname)); 8399566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(resultname,&matmat->resultname)); 8409566063dSJacob Faibussowitsch PetscCall(PetscInfo(A,"Composing %s for product type %s with result %s\n",matmat->composedname,MatProductTypes[matmat->ptype],matmat->resultname ? matmat->resultname : "not specified")); 8419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,matmat->composedname,MatProductSetFromOptions_Shell_X)); 842b77ba244SStefano Zampini PetscFunctionReturn(0); 843b77ba244SStefano Zampini } 844b77ba244SStefano Zampini 845b77ba244SStefano Zampini /*@C 846b77ba244SStefano Zampini MatShellSetMatProductOperation - Allows user to set a matrix matrix operation for a shell matrix. 847b77ba244SStefano Zampini 848b77ba244SStefano Zampini Logically Collective on Mat 849b77ba244SStefano Zampini 850b77ba244SStefano Zampini Input Parameters: 851b77ba244SStefano Zampini + A - the shell matrix 852b77ba244SStefano Zampini . ptype - the product type 853b77ba244SStefano Zampini . symbolic - the function for the symbolic phase (can be NULL) 854b77ba244SStefano Zampini . numeric - the function for the numerical phase 855b77ba244SStefano Zampini . destroy - the function for the destruction of the needed data generated during the symbolic phase (can be NULL) 856b77ba244SStefano Zampini . Btype - the matrix type for the matrix to be multiplied against 857b77ba244SStefano Zampini - Ctype - the matrix type for the result (can be NULL) 858b77ba244SStefano Zampini 859b77ba244SStefano Zampini Level: advanced 860b77ba244SStefano Zampini 861b77ba244SStefano Zampini Usage: 862b77ba244SStefano Zampini $ extern PetscErrorCode usersymbolic(Mat,Mat,Mat,void**); 863b77ba244SStefano Zampini $ extern PetscErrorCode usernumeric(Mat,Mat,Mat,void*); 864b77ba244SStefano Zampini $ extern PetscErrorCode userdestroy(void*); 865b77ba244SStefano Zampini $ MatCreateShell(comm,m,n,M,N,ctx,&A); 866b77ba244SStefano Zampini $ MatShellSetMatProductOperation(A,MATPRODUCT_AB,usersymbolic,usernumeric,userdestroy,MATSEQAIJ,MATDENSE); 867b77ba244SStefano Zampini $ [ create B of type SEQAIJ etc..] 868b77ba244SStefano Zampini $ MatProductCreate(A,B,NULL,&C); 869b77ba244SStefano Zampini $ MatProductSetType(C,MATPRODUCT_AB); 870b77ba244SStefano Zampini $ MatProductSetFromOptions(C); 871b77ba244SStefano Zampini $ MatProductSymbolic(C); -> actually runs the user defined symbolic operation 872b77ba244SStefano Zampini $ MatProductNumeric(C); -> actually runs the user defined numeric operation 873b77ba244SStefano Zampini $ [ use C = A*B ] 874b77ba244SStefano Zampini 875b77ba244SStefano Zampini Notes: 876b77ba244SStefano Zampini MATPRODUCT_ABC is not supported yet. Not supported in Fortran. 877b77ba244SStefano 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. 878b77ba244SStefano Zampini Any additional data needed by the matrix product needs to be returned during the symbolic phase and destroyed with the destroy callback. 879b77ba244SStefano Zampini PETSc will take care of calling the user-defined callbacks. 880b77ba244SStefano Zampini It is allowed to specify the same callbacks for different Btype matrix types. 881b77ba244SStefano Zampini The couple (Btype,ptype) uniquely identifies the operation: the last specified callbacks takes precedence. 882b77ba244SStefano Zampini 883db781477SPatrick Sanan .seealso: `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellSetContext()`, `MatSetOperation()`, `MatProductType`, `MatType`, `MatSetUp()` 884b77ba244SStefano Zampini @*/ 885b77ba244SStefano 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) 886b77ba244SStefano Zampini { 887b77ba244SStefano Zampini PetscFunctionBegin; 888b77ba244SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 889b77ba244SStefano Zampini PetscValidLogicalCollectiveEnum(A,ptype,2); 89008401ef6SPierre Jolivet PetscCheck(ptype != MATPRODUCT_ABC,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not for product type %s",MatProductTypes[ptype]); 89128b400f6SJacob Faibussowitsch PetscCheck(numeric,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Missing numeric routine, argument 4"); 892dadcf809SJacob Faibussowitsch PetscValidCharPointer(Btype,6); 893dadcf809SJacob Faibussowitsch if (Ctype) PetscValidCharPointer(Ctype,7); 894cac4c232SBarry Smith PetscTryMethod(A,"MatShellSetMatProductOperation_C",(Mat,MatProductType,PetscErrorCode(*)(Mat,Mat,Mat,void**),PetscErrorCode(*)(Mat,Mat,Mat,void*),PetscErrorCode(*)(void*),MatType,MatType),(A,ptype,symbolic,numeric,destroy,Btype,Ctype)); 895b77ba244SStefano Zampini PetscFunctionReturn(0); 896b77ba244SStefano Zampini } 897b77ba244SStefano Zampini 898*800f99ffSJeremy L Thompson static PetscErrorCode MatShellSetMatProductOperation_Shell(Mat A,MatProductType ptype,PetscErrorCode (*symbolic)(Mat,Mat,Mat,void**),PetscErrorCode (*numeric)(Mat,Mat,Mat,void*),PetscErrorCode (*destroy)(void *),MatType Btype,MatType Ctype) 899b77ba244SStefano Zampini { 900b77ba244SStefano Zampini PetscBool flg; 901b77ba244SStefano Zampini char composedname[256]; 902b77ba244SStefano Zampini MatRootName Bnames = MatRootNameList, Cnames = MatRootNameList; 903b77ba244SStefano Zampini PetscMPIInt size; 904b77ba244SStefano Zampini 905b77ba244SStefano Zampini PetscFunctionBegin; 906b77ba244SStefano Zampini PetscValidType(A,1); 907b77ba244SStefano Zampini while (Bnames) { /* user passed in the root name */ 9089566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(Btype,Bnames->rname,&flg)); 909b77ba244SStefano Zampini if (flg) break; 910b77ba244SStefano Zampini Bnames = Bnames->next; 911b77ba244SStefano Zampini } 912b77ba244SStefano Zampini while (Cnames) { /* user passed in the root name */ 9139566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(Ctype,Cnames->rname,&flg)); 914b77ba244SStefano Zampini if (flg) break; 915b77ba244SStefano Zampini Cnames = Cnames->next; 916b77ba244SStefano Zampini } 9179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 918b77ba244SStefano Zampini Btype = Bnames ? (size > 1 ? Bnames->mname : Bnames->sname) : Btype; 919b77ba244SStefano Zampini Ctype = Cnames ? (size > 1 ? Cnames->mname : Cnames->sname) : Ctype; 9209566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(composedname,sizeof(composedname),"MatProductSetFromOptions_%s_%s_C",((PetscObject)A)->type_name,Btype)); 9219566063dSJacob Faibussowitsch PetscCall(MatShellSetMatProductOperation_Private(A,ptype,symbolic,numeric,destroy,composedname,Ctype)); 9223a40ed3dSBarry Smith PetscFunctionReturn(0); 923e51e0e81SBarry Smith } 924e51e0e81SBarry Smith 925*800f99ffSJeremy L Thompson static PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str) 9267fabda3fSAlex Fikl { 92728f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data; 928cb8c8a77SBarry Smith PetscBool matflg; 929b77ba244SStefano Zampini MatShellMatFunctionList matmatA; 9307fabda3fSAlex Fikl 9317fabda3fSAlex Fikl PetscFunctionBegin; 9329566063dSJacob Faibussowitsch PetscCall(MatIsShell(B,&matflg)); 93328b400f6SJacob Faibussowitsch PetscCheck(matflg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix %s not derived from MATSHELL",((PetscObject)B)->type_name); 9347fabda3fSAlex Fikl 9359566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps))); 9369566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps))); 9377fabda3fSAlex Fikl 9381baa6e33SBarry Smith if (shellA->ops->copy) PetscCall((*shellA->ops->copy)(A,B,str)); 9397fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 9407fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 9410c0fd78eSBarry Smith if (shellA->dshift) { 9420c0fd78eSBarry Smith if (!shellB->dshift) { 9439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shellA->dshift,&shellB->dshift)); 9447fabda3fSAlex Fikl } 9459566063dSJacob Faibussowitsch PetscCall(VecCopy(shellA->dshift,shellB->dshift)); 9467fabda3fSAlex Fikl } else { 9479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shellB->dshift)); 9487fabda3fSAlex Fikl } 9490c0fd78eSBarry Smith if (shellA->left) { 9500c0fd78eSBarry Smith if (!shellB->left) { 9519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shellA->left,&shellB->left)); 9527fabda3fSAlex Fikl } 9539566063dSJacob Faibussowitsch PetscCall(VecCopy(shellA->left,shellB->left)); 9547fabda3fSAlex Fikl } else { 9559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shellB->left)); 9567fabda3fSAlex Fikl } 9570c0fd78eSBarry Smith if (shellA->right) { 9580c0fd78eSBarry Smith if (!shellB->right) { 9599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shellA->right,&shellB->right)); 9607fabda3fSAlex Fikl } 9619566063dSJacob Faibussowitsch PetscCall(VecCopy(shellA->right,shellB->right)); 9627fabda3fSAlex Fikl } else { 9639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shellB->right)); 9647fabda3fSAlex Fikl } 9659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&shellB->axpy)); 966ef5c7bd2SStefano Zampini shellB->axpy_vscale = 0.0; 967ef5c7bd2SStefano Zampini shellB->axpy_state = 0; 96840e381c3SBarry Smith if (shellA->axpy) { 9699566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)shellA->axpy)); 97040e381c3SBarry Smith shellB->axpy = shellA->axpy; 97140e381c3SBarry Smith shellB->axpy_vscale = shellA->axpy_vscale; 972ef5c7bd2SStefano Zampini shellB->axpy_state = shellA->axpy_state; 97340e381c3SBarry Smith } 97445960306SStefano Zampini if (shellA->zrows) { 9759566063dSJacob Faibussowitsch PetscCall(ISDuplicate(shellA->zrows,&shellB->zrows)); 97645960306SStefano Zampini if (shellA->zcols) { 9779566063dSJacob Faibussowitsch PetscCall(ISDuplicate(shellA->zcols,&shellB->zcols)); 97845960306SStefano Zampini } 9799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shellA->zvals,&shellB->zvals)); 9809566063dSJacob Faibussowitsch PetscCall(VecCopy(shellA->zvals,shellB->zvals)); 9819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shellA->zvals_w,&shellB->zvals_w)); 9829566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)shellA->zvals_sct_r)); 9839566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)shellA->zvals_sct_c)); 98445960306SStefano Zampini shellB->zvals_sct_r = shellA->zvals_sct_r; 98545960306SStefano Zampini shellB->zvals_sct_c = shellA->zvals_sct_c; 98645960306SStefano Zampini } 987b77ba244SStefano Zampini 988b77ba244SStefano Zampini matmatA = shellA->matmat; 989b77ba244SStefano Zampini if (matmatA) { 990b77ba244SStefano Zampini while (matmatA->next) { 9919566063dSJacob Faibussowitsch PetscCall(MatShellSetMatProductOperation_Private(B,matmatA->ptype,matmatA->symbolic,matmatA->numeric,matmatA->destroy,matmatA->composedname,matmatA->resultname)); 992b77ba244SStefano Zampini matmatA = matmatA->next; 993b77ba244SStefano Zampini } 994b77ba244SStefano Zampini } 9957fabda3fSAlex Fikl PetscFunctionReturn(0); 9967fabda3fSAlex Fikl } 9977fabda3fSAlex Fikl 998*800f99ffSJeremy L Thompson static PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M) 999cb8c8a77SBarry Smith { 1000cb8c8a77SBarry Smith PetscFunctionBegin; 1001*800f99ffSJeremy L Thompson PetscCall(MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,NULL,M)); 1002*800f99ffSJeremy L Thompson ((Mat_Shell*)(*M)->data)->ctxcontainer = ((Mat_Shell*)mat->data)->ctxcontainer; 1003*800f99ffSJeremy L Thompson PetscCall(PetscObjectCompose((PetscObject)(*M),"MatShell ctx",(PetscObject)((Mat_Shell*)(*M)->data)->ctxcontainer)); 10049566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)(*M),((PetscObject)mat)->type_name)); 1005a4b1380bSStefano Zampini if (op != MAT_DO_NOT_COPY_VALUES) { 10069566063dSJacob Faibussowitsch PetscCall(MatCopy(mat,*M,SAME_NONZERO_PATTERN)); 1007a4b1380bSStefano Zampini } 1008cb8c8a77SBarry Smith PetscFunctionReturn(0); 1009cb8c8a77SBarry Smith } 1010cb8c8a77SBarry Smith 1011dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 1012ef66eb69SBarry Smith { 1013ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 101425578ef6SJed Brown Vec xx; 1015e3f487b0SBarry Smith PetscObjectState instate,outstate; 1016ef66eb69SBarry Smith 1017ef66eb69SBarry Smith PetscFunctionBegin; 101828b400f6SJacob Faibussowitsch PetscCheck(shell->ops->mult,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL"); 10199566063dSJacob Faibussowitsch PetscCall(MatShellPreZeroRight(A,x,&xx)); 10209566063dSJacob Faibussowitsch PetscCall(MatShellPreScaleRight(A,xx,&xx)); 10219566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 10229566063dSJacob Faibussowitsch PetscCall((*shell->ops->mult)(A,xx,y)); 10239566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 1024e3f487b0SBarry Smith if (instate == outstate) { 1025e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 10269566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)y)); 1027e3f487b0SBarry Smith } 10289566063dSJacob Faibussowitsch PetscCall(MatShellShiftAndScale(A,xx,y)); 10299566063dSJacob Faibussowitsch PetscCall(MatShellPostScaleLeft(A,y)); 10309566063dSJacob Faibussowitsch PetscCall(MatShellPostZeroLeft(A,y)); 10319f137db4SBarry Smith 10329f137db4SBarry Smith if (shell->axpy) { 1033ef5c7bd2SStefano Zampini Mat X; 1034ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 1035ef5c7bd2SStefano Zampini 10369566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(shell->axpy,&X)); 10379566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X,&axpy_state)); 103808401ef6SPierre Jolivet PetscCheck(shell->axpy_state == axpy_state,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)"); 1039b77ba244SStefano Zampini 10409566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(shell->axpy,shell->axpy_right ? NULL : &shell->axpy_right,shell->axpy_left ? NULL : &shell->axpy_left)); 10419566063dSJacob Faibussowitsch PetscCall(VecCopy(x,shell->axpy_right)); 10429566063dSJacob Faibussowitsch PetscCall(MatMult(shell->axpy,shell->axpy_right,shell->axpy_left)); 10439566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,shell->axpy_vscale,shell->axpy_left)); 10449f137db4SBarry Smith } 1045ef66eb69SBarry Smith PetscFunctionReturn(0); 1046ef66eb69SBarry Smith } 1047ef66eb69SBarry Smith 1048*800f99ffSJeremy L Thompson static PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 10495edf6516SJed Brown { 10505edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 10515edf6516SJed Brown 10525edf6516SJed Brown PetscFunctionBegin; 10535edf6516SJed Brown if (y == z) { 10549566063dSJacob Faibussowitsch if (!shell->right_add_work) PetscCall(VecDuplicate(z,&shell->right_add_work)); 10559566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,shell->right_add_work)); 10569566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,shell->right_add_work)); 10575edf6516SJed Brown } else { 10589566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,z)); 10599566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,y)); 10605edf6516SJed Brown } 10615edf6516SJed Brown PetscFunctionReturn(0); 10625edf6516SJed Brown } 10635edf6516SJed Brown 1064*800f99ffSJeremy L Thompson static PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 106518be62a5SSatish Balay { 106618be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 106725578ef6SJed Brown Vec xx; 1068e3f487b0SBarry Smith PetscObjectState instate,outstate; 106918be62a5SSatish Balay 107018be62a5SSatish Balay PetscFunctionBegin; 107128b400f6SJacob Faibussowitsch PetscCheck(shell->ops->multtranspose,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL"); 10729566063dSJacob Faibussowitsch PetscCall(MatShellPreZeroLeft(A,x,&xx)); 10739566063dSJacob Faibussowitsch PetscCall(MatShellPreScaleLeft(A,xx,&xx)); 10749566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &instate)); 10759566063dSJacob Faibussowitsch PetscCall((*shell->ops->multtranspose)(A,xx,y)); 10769566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)y, &outstate)); 1077e3f487b0SBarry Smith if (instate == outstate) { 1078e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 10799566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)y)); 1080e3f487b0SBarry Smith } 10819566063dSJacob Faibussowitsch PetscCall(MatShellShiftAndScale(A,xx,y)); 10829566063dSJacob Faibussowitsch PetscCall(MatShellPostScaleRight(A,y)); 10839566063dSJacob Faibussowitsch PetscCall(MatShellPostZeroRight(A,y)); 1084350ec83bSStefano Zampini 1085350ec83bSStefano Zampini if (shell->axpy) { 1086ef5c7bd2SStefano Zampini Mat X; 1087ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 1088ef5c7bd2SStefano Zampini 10899566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(shell->axpy,&X)); 10909566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X,&axpy_state)); 109108401ef6SPierre Jolivet PetscCheck(shell->axpy_state == axpy_state,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)"); 10929566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(shell->axpy,shell->axpy_right ? NULL : &shell->axpy_right,shell->axpy_left ? NULL : &shell->axpy_left)); 10939566063dSJacob Faibussowitsch PetscCall(VecCopy(x,shell->axpy_left)); 10949566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(shell->axpy,shell->axpy_left,shell->axpy_right)); 10959566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,shell->axpy_vscale,shell->axpy_right)); 1096350ec83bSStefano Zampini } 109718be62a5SSatish Balay PetscFunctionReturn(0); 109818be62a5SSatish Balay } 109918be62a5SSatish Balay 1100*800f99ffSJeremy L Thompson static PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 11015edf6516SJed Brown { 11025edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 11035edf6516SJed Brown 11045edf6516SJed Brown PetscFunctionBegin; 11055edf6516SJed Brown if (y == z) { 11069566063dSJacob Faibussowitsch if (!shell->left_add_work) PetscCall(VecDuplicate(z,&shell->left_add_work)); 11079566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,shell->left_add_work)); 11089566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,shell->left_add_work)); 11095edf6516SJed Brown } else { 11109566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,z)); 11119566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,1.0,y)); 11125edf6516SJed Brown } 11135edf6516SJed Brown PetscFunctionReturn(0); 11145edf6516SJed Brown } 11155edf6516SJed Brown 11160c0fd78eSBarry Smith /* 11170c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 11180c0fd78eSBarry Smith */ 1119*800f99ffSJeremy L Thompson static PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 112018be62a5SSatish Balay { 112118be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 112218be62a5SSatish Balay 112318be62a5SSatish Balay PetscFunctionBegin; 11240c0fd78eSBarry Smith if (shell->ops->getdiagonal) { 11259566063dSJacob Faibussowitsch PetscCall((*shell->ops->getdiagonal)(A,v)); 1126305211d5SBarry 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,...)"); 11279566063dSJacob Faibussowitsch PetscCall(VecScale(v,shell->vscale)); 11281baa6e33SBarry Smith if (shell->dshift) PetscCall(VecAXPY(v,1.0,shell->dshift)); 11299566063dSJacob Faibussowitsch PetscCall(VecShift(v,shell->vshift)); 11309566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseMult(v,v,shell->left)); 11319566063dSJacob Faibussowitsch if (shell->right) PetscCall(VecPointwiseMult(v,v,shell->right)); 113245960306SStefano Zampini if (shell->zrows) { 11339566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE)); 11349566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE)); 113545960306SStefano Zampini } 113681c02519SBarry Smith if (shell->axpy) { 1137ef5c7bd2SStefano Zampini Mat X; 1138ef5c7bd2SStefano Zampini PetscObjectState axpy_state; 1139ef5c7bd2SStefano Zampini 11409566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(shell->axpy,&X)); 11419566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X,&axpy_state)); 114208401ef6SPierre Jolivet PetscCheck(shell->axpy_state == axpy_state,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Invalid AXPY state: cannot modify the X matrix passed to MatAXPY(Y,a,X,...)"); 11439566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(shell->axpy,NULL,shell->axpy_left ? NULL : &shell->axpy_left)); 11449566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(shell->axpy,shell->axpy_left)); 11459566063dSJacob Faibussowitsch PetscCall(VecAXPY(v,shell->axpy_vscale,shell->axpy_left)); 114681c02519SBarry Smith } 114718be62a5SSatish Balay PetscFunctionReturn(0); 114818be62a5SSatish Balay } 114918be62a5SSatish Balay 1150*800f99ffSJeremy L Thompson static PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 1151ef66eb69SBarry Smith { 1152ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 1153789d8953SBarry Smith PetscBool flg; 1154b24ad042SBarry Smith 1155ef66eb69SBarry Smith PetscFunctionBegin; 11569566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(Y,&flg)); 115728b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot shift shell matrix if it is not congruent"); 11580c0fd78eSBarry Smith if (shell->left || shell->right) { 11598fe8eb27SJed Brown if (!shell->dshift) { 11609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift)); 11619566063dSJacob Faibussowitsch PetscCall(VecSet(shell->dshift,a)); 11620c0fd78eSBarry Smith } else { 11639566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseMult(shell->dshift,shell->dshift,shell->left)); 11649566063dSJacob Faibussowitsch if (shell->right) PetscCall(VecPointwiseMult(shell->dshift,shell->dshift,shell->right)); 11659566063dSJacob Faibussowitsch PetscCall(VecShift(shell->dshift,a)); 11660c0fd78eSBarry Smith } 11679566063dSJacob Faibussowitsch if (shell->left) PetscCall(VecPointwiseDivide(shell->dshift,shell->dshift,shell->left)); 11689566063dSJacob Faibussowitsch if (shell->right) PetscCall(VecPointwiseDivide(shell->dshift,shell->dshift,shell->right)); 11698fe8eb27SJed Brown } else shell->vshift += a; 11701baa6e33SBarry Smith if (shell->zrows) PetscCall(VecShift(shell->zvals,a)); 1171ef66eb69SBarry Smith PetscFunctionReturn(0); 1172ef66eb69SBarry Smith } 1173ef66eb69SBarry Smith 1174*800f99ffSJeremy L Thompson static PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s) 11756464bf51SAlex Fikl { 11766464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 11776464bf51SAlex Fikl 11786464bf51SAlex Fikl PetscFunctionBegin; 11799566063dSJacob Faibussowitsch if (!shell->dshift) PetscCall(VecDuplicate(D,&shell->dshift)); 11800c0fd78eSBarry Smith if (shell->left || shell->right) { 11819566063dSJacob Faibussowitsch if (!shell->right_work) PetscCall(VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work)); 11829f137db4SBarry Smith if (shell->left && shell->right) { 11839566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(shell->right_work,D,shell->left)); 11849566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(shell->right_work,shell->right_work,shell->right)); 11859f137db4SBarry Smith } else if (shell->left) { 11869566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(shell->right_work,D,shell->left)); 11879f137db4SBarry Smith } else { 11889566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(shell->right_work,D,shell->right)); 11899f137db4SBarry Smith } 11909566063dSJacob Faibussowitsch PetscCall(VecAXPY(shell->dshift,s,shell->right_work)); 11910c0fd78eSBarry Smith } else { 11929566063dSJacob Faibussowitsch PetscCall(VecAXPY(shell->dshift,s,D)); 1193b253ae0bS“Marek } 1194b253ae0bS“Marek PetscFunctionReturn(0); 1195b253ae0bS“Marek } 1196b253ae0bS“Marek 1197b253ae0bS“Marek PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 1198b253ae0bS“Marek { 119945960306SStefano Zampini Mat_Shell *shell = (Mat_Shell*)A->data; 1200b253ae0bS“Marek Vec d; 1201789d8953SBarry Smith PetscBool flg; 1202b253ae0bS“Marek 1203b253ae0bS“Marek PetscFunctionBegin; 12049566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A,&flg)); 120528b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot diagonal set or shift shell matrix if it is not congruent"); 1206b253ae0bS“Marek if (ins == INSERT_VALUES) { 120728b400f6SJacob Faibussowitsch PetscCheck(A->ops->getdiagonal,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first"); 12089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(D,&d)); 12099566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,d)); 12109566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet_Shell_Private(A,d,-1.)); 12119566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet_Shell_Private(A,D,1.)); 12129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d)); 12131baa6e33SBarry Smith if (shell->zrows) PetscCall(VecCopy(D,shell->zvals)); 1214b253ae0bS“Marek } else { 12159566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet_Shell_Private(A,D,1.)); 12161baa6e33SBarry Smith if (shell->zrows) PetscCall(VecAXPY(shell->zvals,1.0,D)); 12176464bf51SAlex Fikl } 12186464bf51SAlex Fikl PetscFunctionReturn(0); 12196464bf51SAlex Fikl } 12206464bf51SAlex Fikl 1221*800f99ffSJeremy L Thompson static PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 1222ef66eb69SBarry Smith { 1223ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 1224b24ad042SBarry Smith 1225ef66eb69SBarry Smith PetscFunctionBegin; 1226f4df32b1SMatthew Knepley shell->vscale *= a; 12270c0fd78eSBarry Smith shell->vshift *= a; 12281baa6e33SBarry Smith if (shell->dshift) PetscCall(VecScale(shell->dshift,a)); 122981c02519SBarry Smith shell->axpy_vscale *= a; 12301baa6e33SBarry Smith if (shell->zrows) PetscCall(VecScale(shell->zvals,a)); 12318fe8eb27SJed Brown PetscFunctionReturn(0); 123218be62a5SSatish Balay } 12338fe8eb27SJed Brown 12348fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 12358fe8eb27SJed Brown { 12368fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 12378fe8eb27SJed Brown 12388fe8eb27SJed Brown PetscFunctionBegin; 12398fe8eb27SJed Brown if (left) { 12400c0fd78eSBarry Smith if (!shell->left) { 12419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&shell->left)); 12429566063dSJacob Faibussowitsch PetscCall(VecCopy(left,shell->left)); 12430c0fd78eSBarry Smith } else { 12449566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->left,shell->left,left)); 124518be62a5SSatish Balay } 12461baa6e33SBarry Smith if (shell->zrows) PetscCall(VecPointwiseMult(shell->zvals,shell->zvals,left)); 1247ef66eb69SBarry Smith } 12488fe8eb27SJed Brown if (right) { 12490c0fd78eSBarry Smith if (!shell->right) { 12509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right,&shell->right)); 12519566063dSJacob Faibussowitsch PetscCall(VecCopy(right,shell->right)); 12520c0fd78eSBarry Smith } else { 12539566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->right,shell->right,right)); 12548fe8eb27SJed Brown } 125545960306SStefano Zampini if (shell->zrows) { 125645960306SStefano Zampini if (!shell->left_work) { 12579566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Y,NULL,&shell->left_work)); 125845960306SStefano Zampini } 12599566063dSJacob Faibussowitsch PetscCall(VecSet(shell->zvals_w,1.0)); 12609566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 12619566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD)); 12629566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(shell->zvals,shell->zvals,shell->zvals_w)); 126345960306SStefano Zampini } 12648fe8eb27SJed Brown } 12651baa6e33SBarry Smith if (shell->axpy) PetscCall(MatDiagonalScale(shell->axpy,left,right)); 1266ef66eb69SBarry Smith PetscFunctionReturn(0); 1267ef66eb69SBarry Smith } 1268ef66eb69SBarry Smith 1269dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 1270ef66eb69SBarry Smith { 1271ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 1272ef66eb69SBarry Smith 1273ef66eb69SBarry Smith PetscFunctionBegin; 12748fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 1275ef66eb69SBarry Smith shell->vshift = 0.0; 1276ef66eb69SBarry Smith shell->vscale = 1.0; 1277ef5c7bd2SStefano Zampini shell->axpy_vscale = 0.0; 1278ef5c7bd2SStefano Zampini shell->axpy_state = 0; 12799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->dshift)); 12809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->left)); 12819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->right)); 12829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&shell->axpy)); 12839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->axpy_left)); 12849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&shell->axpy_right)); 12859566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_c)); 12869566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&shell->zvals_sct_r)); 12879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zrows)); 12889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&shell->zcols)); 1289ef66eb69SBarry Smith } 1290ef66eb69SBarry Smith PetscFunctionReturn(0); 1291ef66eb69SBarry Smith } 1292ef66eb69SBarry Smith 12933b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 12943b49f96aSBarry Smith { 12953b49f96aSBarry Smith PetscFunctionBegin; 12963b49f96aSBarry Smith *missing = PETSC_FALSE; 12973b49f96aSBarry Smith PetscFunctionReturn(0); 12983b49f96aSBarry Smith } 12993b49f96aSBarry Smith 13009f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str) 13019f137db4SBarry Smith { 13029f137db4SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 13039f137db4SBarry Smith 13049f137db4SBarry Smith PetscFunctionBegin; 1305ef5c7bd2SStefano Zampini if (X == Y) { 13069566063dSJacob Faibussowitsch PetscCall(MatScale(Y,1.0 + a)); 1307ef5c7bd2SStefano Zampini PetscFunctionReturn(0); 1308ef5c7bd2SStefano Zampini } 1309ef5c7bd2SStefano Zampini if (!shell->axpy) { 13109566063dSJacob Faibussowitsch PetscCall(MatConvertFrom_Shell(X,MATSHELL,MAT_INITIAL_MATRIX,&shell->axpy)); 13119f137db4SBarry Smith shell->axpy_vscale = a; 13129566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)X,&shell->axpy_state)); 1313ef5c7bd2SStefano Zampini } else { 13149566063dSJacob Faibussowitsch PetscCall(MatAXPY(shell->axpy,a/shell->axpy_vscale,X,str)); 1315ef5c7bd2SStefano Zampini } 13169f137db4SBarry Smith PetscFunctionReturn(0); 13179f137db4SBarry Smith } 13189f137db4SBarry Smith 1319f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL, 1320f4259b30SLisandro Dalcin NULL, 1321f4259b30SLisandro Dalcin NULL, 1322f4259b30SLisandro Dalcin NULL, 13230c0fd78eSBarry Smith /* 4*/ MatMultAdd_Shell, 1324f4259b30SLisandro Dalcin NULL, 13250c0fd78eSBarry Smith MatMultTransposeAdd_Shell, 1326f4259b30SLisandro Dalcin NULL, 1327f4259b30SLisandro Dalcin NULL, 1328f4259b30SLisandro Dalcin NULL, 1329f4259b30SLisandro Dalcin /*10*/ NULL, 1330f4259b30SLisandro Dalcin NULL, 1331f4259b30SLisandro Dalcin NULL, 1332f4259b30SLisandro Dalcin NULL, 1333f4259b30SLisandro Dalcin NULL, 1334f4259b30SLisandro Dalcin /*15*/ NULL, 1335f4259b30SLisandro Dalcin NULL, 1336f4259b30SLisandro Dalcin NULL, 13378fe8eb27SJed Brown MatDiagonalScale_Shell, 1338f4259b30SLisandro Dalcin NULL, 1339f4259b30SLisandro Dalcin /*20*/ NULL, 1340ef66eb69SBarry Smith MatAssemblyEnd_Shell, 1341f4259b30SLisandro Dalcin NULL, 1342f4259b30SLisandro Dalcin NULL, 134345960306SStefano Zampini /*24*/ MatZeroRows_Shell, 1344f4259b30SLisandro Dalcin NULL, 1345f4259b30SLisandro Dalcin NULL, 1346f4259b30SLisandro Dalcin NULL, 1347f4259b30SLisandro Dalcin NULL, 1348f4259b30SLisandro Dalcin /*29*/ NULL, 1349f4259b30SLisandro Dalcin NULL, 1350f4259b30SLisandro Dalcin NULL, 1351f4259b30SLisandro Dalcin NULL, 1352f4259b30SLisandro Dalcin NULL, 1353cb8c8a77SBarry Smith /*34*/ MatDuplicate_Shell, 1354f4259b30SLisandro Dalcin NULL, 1355f4259b30SLisandro Dalcin NULL, 1356f4259b30SLisandro Dalcin NULL, 1357f4259b30SLisandro Dalcin NULL, 13589f137db4SBarry Smith /*39*/ MatAXPY_Shell, 1359f4259b30SLisandro Dalcin NULL, 1360f4259b30SLisandro Dalcin NULL, 1361f4259b30SLisandro Dalcin NULL, 1362cb8c8a77SBarry Smith MatCopy_Shell, 1363f4259b30SLisandro Dalcin /*44*/ NULL, 1364ef66eb69SBarry Smith MatScale_Shell, 1365ef66eb69SBarry Smith MatShift_Shell, 13666464bf51SAlex Fikl MatDiagonalSet_Shell, 136745960306SStefano Zampini MatZeroRowsColumns_Shell, 1368f4259b30SLisandro Dalcin /*49*/ NULL, 1369f4259b30SLisandro Dalcin NULL, 1370f4259b30SLisandro Dalcin NULL, 1371f4259b30SLisandro Dalcin NULL, 1372f4259b30SLisandro Dalcin NULL, 1373f4259b30SLisandro Dalcin /*54*/ NULL, 1374f4259b30SLisandro Dalcin NULL, 1375f4259b30SLisandro Dalcin NULL, 1376f4259b30SLisandro Dalcin NULL, 1377f4259b30SLisandro Dalcin NULL, 1378f4259b30SLisandro Dalcin /*59*/ NULL, 1379b9b97703SBarry Smith MatDestroy_Shell, 1380f4259b30SLisandro Dalcin NULL, 1381251fad3fSStefano Zampini MatConvertFrom_Shell, 1382f4259b30SLisandro Dalcin NULL, 1383f4259b30SLisandro Dalcin /*64*/ NULL, 1384f4259b30SLisandro Dalcin NULL, 1385f4259b30SLisandro Dalcin NULL, 1386f4259b30SLisandro Dalcin NULL, 1387f4259b30SLisandro Dalcin NULL, 1388f4259b30SLisandro Dalcin /*69*/ NULL, 1389f4259b30SLisandro Dalcin NULL, 1390c87e5d42SMatthew Knepley MatConvert_Shell, 1391f4259b30SLisandro Dalcin NULL, 1392f4259b30SLisandro Dalcin NULL, 1393f4259b30SLisandro Dalcin /*74*/ NULL, 1394f4259b30SLisandro Dalcin NULL, 1395f4259b30SLisandro Dalcin NULL, 1396f4259b30SLisandro Dalcin NULL, 1397f4259b30SLisandro Dalcin NULL, 1398f4259b30SLisandro Dalcin /*79*/ NULL, 1399f4259b30SLisandro Dalcin NULL, 1400f4259b30SLisandro Dalcin NULL, 1401f4259b30SLisandro Dalcin NULL, 1402f4259b30SLisandro Dalcin NULL, 1403f4259b30SLisandro Dalcin /*84*/ NULL, 1404f4259b30SLisandro Dalcin NULL, 1405f4259b30SLisandro Dalcin NULL, 1406f4259b30SLisandro Dalcin NULL, 1407f4259b30SLisandro Dalcin NULL, 1408f4259b30SLisandro Dalcin /*89*/ NULL, 1409f4259b30SLisandro Dalcin NULL, 1410f4259b30SLisandro Dalcin NULL, 1411f4259b30SLisandro Dalcin NULL, 1412f4259b30SLisandro Dalcin NULL, 1413f4259b30SLisandro Dalcin /*94*/ NULL, 1414f4259b30SLisandro Dalcin NULL, 1415f4259b30SLisandro Dalcin NULL, 1416f4259b30SLisandro Dalcin NULL, 1417f4259b30SLisandro Dalcin NULL, 1418f4259b30SLisandro Dalcin /*99*/ NULL, 1419f4259b30SLisandro Dalcin NULL, 1420f4259b30SLisandro Dalcin NULL, 1421f4259b30SLisandro Dalcin NULL, 1422f4259b30SLisandro Dalcin NULL, 1423f4259b30SLisandro Dalcin /*104*/ NULL, 1424f4259b30SLisandro Dalcin NULL, 1425f4259b30SLisandro Dalcin NULL, 1426f4259b30SLisandro Dalcin NULL, 1427f4259b30SLisandro Dalcin NULL, 1428f4259b30SLisandro Dalcin /*109*/ NULL, 1429f4259b30SLisandro Dalcin NULL, 1430f4259b30SLisandro Dalcin NULL, 1431f4259b30SLisandro Dalcin NULL, 14323b49f96aSBarry Smith MatMissingDiagonal_Shell, 1433f4259b30SLisandro Dalcin /*114*/ NULL, 1434f4259b30SLisandro Dalcin NULL, 1435f4259b30SLisandro Dalcin NULL, 1436f4259b30SLisandro Dalcin NULL, 1437f4259b30SLisandro Dalcin NULL, 1438f4259b30SLisandro Dalcin /*119*/ NULL, 1439f4259b30SLisandro Dalcin NULL, 1440f4259b30SLisandro Dalcin NULL, 1441f4259b30SLisandro Dalcin NULL, 1442f4259b30SLisandro Dalcin NULL, 1443f4259b30SLisandro Dalcin /*124*/ NULL, 1444f4259b30SLisandro Dalcin NULL, 1445f4259b30SLisandro Dalcin NULL, 1446f4259b30SLisandro Dalcin NULL, 1447f4259b30SLisandro Dalcin NULL, 1448f4259b30SLisandro Dalcin /*129*/ NULL, 1449f4259b30SLisandro Dalcin NULL, 1450f4259b30SLisandro Dalcin NULL, 1451f4259b30SLisandro Dalcin NULL, 1452f4259b30SLisandro Dalcin NULL, 1453f4259b30SLisandro Dalcin /*134*/ NULL, 1454f4259b30SLisandro Dalcin NULL, 1455f4259b30SLisandro Dalcin NULL, 1456f4259b30SLisandro Dalcin NULL, 1457f4259b30SLisandro Dalcin NULL, 1458f4259b30SLisandro Dalcin /*139*/ NULL, 1459f4259b30SLisandro Dalcin NULL, 1460d70f29a3SPierre Jolivet NULL, 1461d70f29a3SPierre Jolivet NULL, 1462d70f29a3SPierre Jolivet NULL, 1463d70f29a3SPierre Jolivet /*144*/ NULL, 1464d70f29a3SPierre Jolivet NULL, 1465d70f29a3SPierre Jolivet NULL, 146699a7f59eSMark Adams NULL, 146799a7f59eSMark Adams NULL, 14687fb60732SBarry Smith NULL, 14697fb60732SBarry Smith /*150*/ NULL 14703964eb88SJed Brown }; 1471273d9f13SBarry Smith 1472*800f99ffSJeremy L Thompson static PetscErrorCode MatShellSetContext_Shell(Mat mat,void *ctx) 1473789d8953SBarry Smith { 1474789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1475789d8953SBarry Smith 1476789d8953SBarry Smith PetscFunctionBegin; 1477*800f99ffSJeremy L Thompson if (ctx) { 1478*800f99ffSJeremy L Thompson PetscContainer ctxcontainer; 1479*800f99ffSJeremy L Thompson PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)mat),&ctxcontainer)); 1480*800f99ffSJeremy L Thompson PetscCall(PetscContainerSetPointer(ctxcontainer,ctx)); 1481*800f99ffSJeremy L Thompson PetscCall(PetscObjectCompose((PetscObject)mat,"MatShell ctx",(PetscObject)ctxcontainer)); 1482*800f99ffSJeremy L Thompson shell->ctxcontainer = ctxcontainer; 1483*800f99ffSJeremy L Thompson PetscCall(PetscContainerDestroy(&ctxcontainer)); 1484*800f99ffSJeremy L Thompson } else { 1485*800f99ffSJeremy L Thompson PetscCall(PetscObjectCompose((PetscObject)mat,"MatShell ctx",NULL)); 1486*800f99ffSJeremy L Thompson shell->ctxcontainer = NULL; 1487*800f99ffSJeremy L Thompson } 1488*800f99ffSJeremy L Thompson 1489*800f99ffSJeremy L Thompson PetscFunctionReturn(0); 1490*800f99ffSJeremy L Thompson } 1491*800f99ffSJeremy L Thompson 1492*800f99ffSJeremy L Thompson PetscErrorCode MatShellSetContextDestroy_Shell(Mat mat,PetscErrorCode (*f)(void*)) 1493*800f99ffSJeremy L Thompson { 1494*800f99ffSJeremy L Thompson Mat_Shell *shell = (Mat_Shell*)mat->data; 1495*800f99ffSJeremy L Thompson 1496*800f99ffSJeremy L Thompson PetscFunctionBegin; 1497*800f99ffSJeremy L Thompson if (shell->ctxcontainer) PetscCall(PetscContainerSetUserDestroy(shell->ctxcontainer,f)); 1498789d8953SBarry Smith PetscFunctionReturn(0); 1499789d8953SBarry Smith } 1500789d8953SBarry Smith 1501db77db73SJeremy L Thompson static PetscErrorCode MatShellSetVecType_Shell(Mat mat,VecType vtype) 1502db77db73SJeremy L Thompson { 1503db77db73SJeremy L Thompson PetscFunctionBegin; 15049566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->defaultvectype)); 15059566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(vtype,(char**)&mat->defaultvectype)); 1506db77db73SJeremy L Thompson PetscFunctionReturn(0); 1507db77db73SJeremy L Thompson } 1508db77db73SJeremy L Thompson 1509789d8953SBarry Smith PetscErrorCode MatShellSetManageScalingShifts_Shell(Mat A) 1510789d8953SBarry Smith { 1511789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 1512789d8953SBarry Smith 1513789d8953SBarry Smith PetscFunctionBegin; 1514789d8953SBarry Smith shell->managescalingshifts = PETSC_FALSE; 1515789d8953SBarry Smith A->ops->diagonalset = NULL; 1516789d8953SBarry Smith A->ops->diagonalscale = NULL; 1517789d8953SBarry Smith A->ops->scale = NULL; 1518789d8953SBarry Smith A->ops->shift = NULL; 1519789d8953SBarry Smith A->ops->axpy = NULL; 1520789d8953SBarry Smith PetscFunctionReturn(0); 1521789d8953SBarry Smith } 1522789d8953SBarry Smith 1523*800f99ffSJeremy L Thompson static PetscErrorCode MatShellSetOperation_Shell(Mat mat,MatOperation op,void (*f)(void)) 1524789d8953SBarry Smith { 1525feb237baSPierre Jolivet Mat_Shell *shell = (Mat_Shell*)mat->data; 1526789d8953SBarry Smith 1527789d8953SBarry Smith PetscFunctionBegin; 1528789d8953SBarry Smith switch (op) { 1529789d8953SBarry Smith case MATOP_DESTROY: 1530789d8953SBarry Smith shell->ops->destroy = (PetscErrorCode (*)(Mat))f; 1531789d8953SBarry Smith break; 1532789d8953SBarry Smith case MATOP_VIEW: 1533789d8953SBarry Smith if (!mat->ops->viewnative) { 1534789d8953SBarry Smith mat->ops->viewnative = mat->ops->view; 1535789d8953SBarry Smith } 1536789d8953SBarry Smith mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 1537789d8953SBarry Smith break; 1538789d8953SBarry Smith case MATOP_COPY: 1539789d8953SBarry Smith shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f; 1540789d8953SBarry Smith break; 1541789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1542789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1543789d8953SBarry Smith case MATOP_SHIFT: 1544789d8953SBarry Smith case MATOP_SCALE: 1545789d8953SBarry Smith case MATOP_AXPY: 1546789d8953SBarry Smith case MATOP_ZERO_ROWS: 1547789d8953SBarry Smith case MATOP_ZERO_ROWS_COLUMNS: 154828b400f6SJacob Faibussowitsch PetscCheck(!shell->managescalingshifts,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()"); 1549789d8953SBarry Smith (((void(**)(void))mat->ops)[op]) = f; 1550789d8953SBarry Smith break; 1551789d8953SBarry Smith case MATOP_GET_DIAGONAL: 1552789d8953SBarry Smith if (shell->managescalingshifts) { 1553789d8953SBarry Smith shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1554789d8953SBarry Smith mat->ops->getdiagonal = MatGetDiagonal_Shell; 1555789d8953SBarry Smith } else { 1556789d8953SBarry Smith shell->ops->getdiagonal = NULL; 1557789d8953SBarry Smith mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1558789d8953SBarry Smith } 1559789d8953SBarry Smith break; 1560789d8953SBarry Smith case MATOP_MULT: 1561789d8953SBarry Smith if (shell->managescalingshifts) { 1562789d8953SBarry Smith shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1563789d8953SBarry Smith mat->ops->mult = MatMult_Shell; 1564789d8953SBarry Smith } else { 1565789d8953SBarry Smith shell->ops->mult = NULL; 1566789d8953SBarry Smith mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1567789d8953SBarry Smith } 1568789d8953SBarry Smith break; 1569789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 1570789d8953SBarry Smith if (shell->managescalingshifts) { 1571789d8953SBarry Smith shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1572789d8953SBarry Smith mat->ops->multtranspose = MatMultTranspose_Shell; 1573789d8953SBarry Smith } else { 1574789d8953SBarry Smith shell->ops->multtranspose = NULL; 1575789d8953SBarry Smith mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1576789d8953SBarry Smith } 1577789d8953SBarry Smith break; 1578789d8953SBarry Smith default: 1579789d8953SBarry Smith (((void(**)(void))mat->ops)[op]) = f; 1580789d8953SBarry Smith break; 1581789d8953SBarry Smith } 1582789d8953SBarry Smith PetscFunctionReturn(0); 1583789d8953SBarry Smith } 1584789d8953SBarry Smith 1585*800f99ffSJeremy L Thompson static PetscErrorCode MatShellGetOperation_Shell(Mat mat,MatOperation op,void(**f)(void)) 1586789d8953SBarry Smith { 1587789d8953SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1588789d8953SBarry Smith 1589789d8953SBarry Smith PetscFunctionBegin; 1590789d8953SBarry Smith switch (op) { 1591789d8953SBarry Smith case MATOP_DESTROY: 1592789d8953SBarry Smith *f = (void (*)(void))shell->ops->destroy; 1593789d8953SBarry Smith break; 1594789d8953SBarry Smith case MATOP_VIEW: 1595789d8953SBarry Smith *f = (void (*)(void))mat->ops->view; 1596789d8953SBarry Smith break; 1597789d8953SBarry Smith case MATOP_COPY: 1598789d8953SBarry Smith *f = (void (*)(void))shell->ops->copy; 1599789d8953SBarry Smith break; 1600789d8953SBarry Smith case MATOP_DIAGONAL_SET: 1601789d8953SBarry Smith case MATOP_DIAGONAL_SCALE: 1602789d8953SBarry Smith case MATOP_SHIFT: 1603789d8953SBarry Smith case MATOP_SCALE: 1604789d8953SBarry Smith case MATOP_AXPY: 1605789d8953SBarry Smith case MATOP_ZERO_ROWS: 1606789d8953SBarry Smith case MATOP_ZERO_ROWS_COLUMNS: 1607789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1608789d8953SBarry Smith break; 1609789d8953SBarry Smith case MATOP_GET_DIAGONAL: 1610789d8953SBarry Smith if (shell->ops->getdiagonal) 1611789d8953SBarry Smith *f = (void (*)(void))shell->ops->getdiagonal; 1612789d8953SBarry Smith else 1613789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1614789d8953SBarry Smith break; 1615789d8953SBarry Smith case MATOP_MULT: 1616789d8953SBarry Smith if (shell->ops->mult) 1617789d8953SBarry Smith *f = (void (*)(void))shell->ops->mult; 1618789d8953SBarry Smith else 1619789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1620789d8953SBarry Smith break; 1621789d8953SBarry Smith case MATOP_MULT_TRANSPOSE: 1622789d8953SBarry Smith if (shell->ops->multtranspose) 1623789d8953SBarry Smith *f = (void (*)(void))shell->ops->multtranspose; 1624789d8953SBarry Smith else 1625789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1626789d8953SBarry Smith break; 1627789d8953SBarry Smith default: 1628789d8953SBarry Smith *f = (((void (**)(void))mat->ops)[op]); 1629789d8953SBarry Smith } 1630789d8953SBarry Smith PetscFunctionReturn(0); 1631789d8953SBarry Smith } 1632789d8953SBarry Smith 16330bad9183SKris Buschelman /*MC 1634fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 16350bad9183SKris Buschelman 16360bad9183SKris Buschelman Level: advanced 16370bad9183SKris Buschelman 1638db781477SPatrick Sanan .seealso: `MatCreateShell()` 16390bad9183SKris Buschelman M*/ 16400bad9183SKris Buschelman 16418cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 1642273d9f13SBarry Smith { 1643273d9f13SBarry Smith Mat_Shell *b; 1644273d9f13SBarry Smith 1645273d9f13SBarry Smith PetscFunctionBegin; 16469566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps))); 1647273d9f13SBarry Smith 16489566063dSJacob Faibussowitsch PetscCall(PetscNewLog(A,&b)); 1649273d9f13SBarry Smith A->data = (void*)b; 1650273d9f13SBarry Smith 1651*800f99ffSJeremy L Thompson b->ctxcontainer = NULL; 1652ef66eb69SBarry Smith b->vshift = 0.0; 1653ef66eb69SBarry Smith b->vscale = 1.0; 16540c0fd78eSBarry Smith b->managescalingshifts = PETSC_TRUE; 1655273d9f13SBarry Smith A->assembled = PETSC_TRUE; 1656210f0121SBarry Smith A->preallocated = PETSC_FALSE; 16572205254eSKarl Rupp 16589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatShellGetContext_C",MatShellGetContext_Shell)); 16599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatShellSetContext_C",MatShellSetContext_Shell)); 1660*800f99ffSJeremy L Thompson PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatShellSetContextDestroy_C",MatShellSetContextDestroy_Shell)); 16619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatShellSetVecType_C",MatShellSetVecType_Shell)); 16629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatShellSetManageScalingShifts_C",MatShellSetManageScalingShifts_Shell)); 16639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatShellSetOperation_C",MatShellSetOperation_Shell)); 16649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatShellGetOperation_C",MatShellGetOperation_Shell)); 16659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatShellSetMatProductOperation_C",MatShellSetMatProductOperation_Shell)); 16669566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATSHELL)); 1667273d9f13SBarry Smith PetscFunctionReturn(0); 1668273d9f13SBarry Smith } 1669e51e0e81SBarry Smith 16704b828684SBarry Smith /*@C 1671052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 1672ff756334SLois Curfman McInnes private data storage format. 1673e51e0e81SBarry Smith 1674d083f849SBarry Smith Collective 1675c7fcc2eaSBarry Smith 1676e51e0e81SBarry Smith Input Parameters: 1677c7fcc2eaSBarry Smith + comm - MPI communicator 1678c7fcc2eaSBarry Smith . m - number of local rows (must be given) 1679c7fcc2eaSBarry Smith . n - number of local columns (must be given) 1680c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 1681c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 1682c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 1683e51e0e81SBarry Smith 1684ff756334SLois Curfman McInnes Output Parameter: 168544cd7ae7SLois Curfman McInnes . A - the matrix 1686e51e0e81SBarry Smith 1687ff2fd236SBarry Smith Level: advanced 1688ff2fd236SBarry Smith 1689f39d1f56SLois Curfman McInnes Usage: 16905bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1691f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 1692c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1693f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 1694f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 1695f39d1f56SLois Curfman McInnes 1696ff756334SLois Curfman McInnes Notes: 1697ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 1698ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 1699ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 1700e51e0e81SBarry Smith 170195452b02SPatrick Sanan Fortran Notes: 170295452b02SPatrick Sanan To use this from Fortran with a ctx you must write an interface definition for this 1703daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 1704daf670e6SBarry Smith in as the ctx argument. 1705f38a8d56SBarry Smith 1706f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 1707f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 1708645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 1709645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 1710645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 1711645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 1712645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 1713645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 1714645985a0SLois Curfman McInnes For example, 1715f39d1f56SLois Curfman McInnes 1716f39d1f56SLois Curfman McInnes $ 1717f39d1f56SLois Curfman McInnes $ Vec x, y 17185bd1e576SStefano Zampini $ extern PetscErrorCode mult(Mat,Vec,Vec); 1719645985a0SLois Curfman McInnes $ Mat A 1720f39d1f56SLois Curfman McInnes $ 1721c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 1722c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 1723f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 1724c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 1725c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 1726c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 1727645985a0SLois Curfman McInnes $ MatMult(A,x,y); 172845960306SStefano Zampini $ MatDestroy(&A); 172945960306SStefano Zampini $ VecDestroy(&y); 173045960306SStefano Zampini $ VecDestroy(&x); 1731645985a0SLois Curfman McInnes $ 1732e51e0e81SBarry Smith 173345960306SStefano Zampini MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), MatScale(), MatZeroRows() and MatZeroRowsColumns() internally, so these 17349b53d762SBarry Smith operations cannot be overwritten unless MatShellSetManageScalingShifts() is called. 17359b53d762SBarry Smith 17360c0fd78eSBarry Smith For rectangular matrices do all the scalings and shifts make sense? 17370c0fd78eSBarry Smith 173895452b02SPatrick Sanan Developers Notes: 173995452b02SPatrick Sanan Regarding shifting and scaling. The general form is 17400c0fd78eSBarry Smith 17410c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 17420c0fd78eSBarry Smith 17430c0fd78eSBarry Smith The order you apply the operations is important. For example if you have a dshift then 17440c0fd78eSBarry Smith apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift 17450c0fd78eSBarry Smith you get s*vscale*A + diag(shift) 17460c0fd78eSBarry Smith 17470c0fd78eSBarry Smith A is the user provided function. 17480c0fd78eSBarry Smith 1749ad2f5c3fSBarry Smith KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for 1750ad2f5c3fSBarry Smith for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger 1751ad2f5c3fSBarry Smith an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat); 1752ad2f5c3fSBarry Smith each time the MATSHELL matrix has changed. 1753ad2f5c3fSBarry Smith 17547301b172SPierre Jolivet Matrix product operations (i.e. MatMat, MatTransposeMat etc) can be specified using MatShellSetMatProductOperation() 1755b77ba244SStefano Zampini 1756ad2f5c3fSBarry Smith Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided 1757ad2f5c3fSBarry Smith with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale(). 1758ad2f5c3fSBarry Smith 1759db781477SPatrick Sanan .seealso: `MatShellSetOperation()`, `MatHasOperation()`, `MatShellGetContext()`, `MatShellSetContext()`, `MATSHELL`, `MatShellSetManageScalingShifts()`, `MatShellSetMatProductOperation()` 1760e51e0e81SBarry Smith @*/ 17617087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 1762e51e0e81SBarry Smith { 17633a40ed3dSBarry Smith PetscFunctionBegin; 17649566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,A)); 17659566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A,m,n,M,N)); 17669566063dSJacob Faibussowitsch PetscCall(MatSetType(*A,MATSHELL)); 17679566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(*A,ctx)); 17689566063dSJacob Faibussowitsch PetscCall(MatSetUp(*A)); 1769273d9f13SBarry Smith PetscFunctionReturn(0); 1770c7fcc2eaSBarry Smith } 1771c7fcc2eaSBarry Smith 1772c6866cfdSSatish Balay /*@ 1773273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 1774c7fcc2eaSBarry Smith 17753f9fe445SBarry Smith Logically Collective on Mat 1776c7fcc2eaSBarry Smith 1777273d9f13SBarry Smith Input Parameters: 1778273d9f13SBarry Smith + mat - the shell matrix 1779273d9f13SBarry Smith - ctx - the context 1780273d9f13SBarry Smith 1781273d9f13SBarry Smith Level: advanced 1782273d9f13SBarry Smith 178395452b02SPatrick Sanan Fortran Notes: 178495452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 1785daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 1786273d9f13SBarry Smith 1787db781477SPatrick Sanan .seealso: `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()` 17880bc0a719SBarry Smith @*/ 17897087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 1790273d9f13SBarry Smith { 1791273d9f13SBarry Smith PetscFunctionBegin; 17920700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1793cac4c232SBarry Smith PetscTryMethod(mat,"MatShellSetContext_C",(Mat,void*),(mat,ctx)); 17943a40ed3dSBarry Smith PetscFunctionReturn(0); 1795e51e0e81SBarry Smith } 1796e51e0e81SBarry Smith 1797*800f99ffSJeremy L Thompson /*@ 1798*800f99ffSJeremy L Thompson MatShellSetContextDestroy - sets the destroy function for a shell matrix context 1799*800f99ffSJeremy L Thompson 1800*800f99ffSJeremy L Thompson Logically Collective on Mat 1801*800f99ffSJeremy L Thompson 1802*800f99ffSJeremy L Thompson Input Parameters: 1803*800f99ffSJeremy L Thompson + mat - the shell matrix 1804*800f99ffSJeremy L Thompson - f - the context destroy function 1805*800f99ffSJeremy L Thompson 1806*800f99ffSJeremy L Thompson Note: 1807*800f99ffSJeremy L Thompson If the `MatShell` is never duplicated, the behavior of this function is equivalent 1808*800f99ffSJeremy L Thompson to `MatShellSetOperation(Mat,MATOP_DESTROY,f)`. However, `MatShellSetContextDestroy()` 1809*800f99ffSJeremy L Thompson ensures proper reference counting for the user provided context data in the case that 1810*800f99ffSJeremy L Thompson the `MatShell` is duplicated. 1811*800f99ffSJeremy L Thompson 1812*800f99ffSJeremy L Thompson Level: advanced 1813*800f99ffSJeremy L Thompson 1814*800f99ffSJeremy L Thompson .seealso: `MatCreateShell()`, `MatShellSetContext()` 1815*800f99ffSJeremy L Thompson @*/ 1816*800f99ffSJeremy L Thompson PetscErrorCode MatShellSetContextDestroy(Mat mat,PetscErrorCode (*f)(void*)) 1817*800f99ffSJeremy L Thompson { 1818*800f99ffSJeremy L Thompson PetscFunctionBegin; 1819*800f99ffSJeremy L Thompson PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1820*800f99ffSJeremy L Thompson PetscTryMethod(mat,"MatShellSetContextDestroy_C",(Mat,PetscErrorCode (*)(void*)),(mat,f)); 1821*800f99ffSJeremy L Thompson PetscFunctionReturn(0); 1822*800f99ffSJeremy L Thompson } 1823*800f99ffSJeremy L Thompson 1824db77db73SJeremy L Thompson /*@C 1825db77db73SJeremy L Thompson MatShellSetVecType - Sets the type of Vec returned by MatCreateVecs() 1826db77db73SJeremy L Thompson 1827db77db73SJeremy L Thompson Logically collective 1828db77db73SJeremy L Thompson 1829db77db73SJeremy L Thompson Input Parameters: 1830db77db73SJeremy L Thompson + mat - the shell matrix 1831db77db73SJeremy L Thompson - vtype - type to use for creating vectors 1832db77db73SJeremy L Thompson 1833db77db73SJeremy L Thompson Notes: 1834db77db73SJeremy L Thompson 1835db77db73SJeremy L Thompson Level: advanced 1836db77db73SJeremy L Thompson 1837db781477SPatrick Sanan .seealso: `MatCreateVecs()` 1838db77db73SJeremy L Thompson @*/ 1839db77db73SJeremy L Thompson PetscErrorCode MatShellSetVecType(Mat mat,VecType vtype) 1840db77db73SJeremy L Thompson { 1841db77db73SJeremy L Thompson PetscFunctionBegin; 1842cac4c232SBarry Smith PetscTryMethod(mat,"MatShellSetVecType_C",(Mat,VecType),(mat,vtype)); 1843db77db73SJeremy L Thompson PetscFunctionReturn(0); 1844db77db73SJeremy L Thompson } 1845db77db73SJeremy L Thompson 18460c0fd78eSBarry Smith /*@ 18470c0fd78eSBarry Smith MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately 18480c0fd78eSBarry Smith after MatCreateShell() 18490c0fd78eSBarry Smith 18500c0fd78eSBarry Smith Logically Collective on Mat 18510c0fd78eSBarry Smith 18520c0fd78eSBarry Smith Input Parameter: 18530c0fd78eSBarry Smith . mat - the shell matrix 18540c0fd78eSBarry Smith 18550c0fd78eSBarry Smith Level: advanced 18560c0fd78eSBarry Smith 1857db781477SPatrick Sanan .seealso: `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellSetContext()`, `MatShellSetOperation()` 18580c0fd78eSBarry Smith @*/ 18590c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A) 18600c0fd78eSBarry Smith { 18610c0fd78eSBarry Smith PetscFunctionBegin; 18620c0fd78eSBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1863cac4c232SBarry Smith PetscTryMethod(A,"MatShellSetManageScalingShifts_C",(Mat),(A)); 18640c0fd78eSBarry Smith PetscFunctionReturn(0); 18650c0fd78eSBarry Smith } 18660c0fd78eSBarry Smith 1867c16cb8f2SBarry Smith /*@C 1868f3b1f45cSBarry Smith MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function. 1869f3b1f45cSBarry Smith 1870f3b1f45cSBarry Smith Logically Collective on Mat 1871f3b1f45cSBarry Smith 1872f3b1f45cSBarry Smith Input Parameters: 1873f3b1f45cSBarry Smith + mat - the shell matrix 1874f3b1f45cSBarry Smith . f - the function 1875f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1876f3b1f45cSBarry Smith - ctx - an optional context for the function 1877f3b1f45cSBarry Smith 1878f3b1f45cSBarry Smith Output Parameter: 1879f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1880f3b1f45cSBarry Smith 1881f3b1f45cSBarry Smith Options Database: 1882f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1883f3b1f45cSBarry Smith 1884f3b1f45cSBarry Smith Level: advanced 1885f3b1f45cSBarry Smith 188695452b02SPatrick Sanan Fortran Notes: 188795452b02SPatrick Sanan Not supported from Fortran 1888f3b1f45cSBarry Smith 1889db781477SPatrick Sanan .seealso: `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellTestMultTranspose()` 1890f3b1f45cSBarry Smith @*/ 1891f3b1f45cSBarry Smith PetscErrorCode MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1892f3b1f45cSBarry Smith { 1893f3b1f45cSBarry Smith PetscInt m,n; 1894f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1895f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 189674e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1897f3b1f45cSBarry Smith 1898f3b1f45cSBarry Smith PetscFunctionBegin; 1899f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 19009566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v)); 19019566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat,&m,&n)); 19029566063dSJacob Faibussowitsch PetscCall(MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf)); 19039566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(mf,f,ctx)); 19049566063dSJacob Faibussowitsch PetscCall(MatMFFDSetBase(mf,base,NULL)); 1905f3b1f45cSBarry Smith 19069566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(mf,MATAIJ,&Dmf)); 19079566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(mat,MATAIJ,&Dmat)); 1908f3b1f45cSBarry Smith 19099566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff)); 19109566063dSJacob Faibussowitsch PetscCall(MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN)); 19119566063dSJacob Faibussowitsch PetscCall(MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm)); 19129566063dSJacob Faibussowitsch PetscCall(MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm)); 1913f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1914f3b1f45cSBarry Smith flag = PETSC_FALSE; 1915f3b1f45cSBarry Smith if (v) { 19169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm))); 19179566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view")); 19189566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view")); 19199566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view")); 1920f3b1f45cSBarry Smith } 1921f3b1f45cSBarry Smith } else if (v) { 19229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n")); 1923f3b1f45cSBarry Smith } 1924f3b1f45cSBarry Smith if (flg) *flg = flag; 19259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ddiff)); 19269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mf)); 19279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Dmf)); 19289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Dmat)); 1929f3b1f45cSBarry Smith PetscFunctionReturn(0); 1930f3b1f45cSBarry Smith } 1931f3b1f45cSBarry Smith 1932f3b1f45cSBarry Smith /*@C 19337301b172SPierre Jolivet MatShellTestMultTranspose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function. 1934f3b1f45cSBarry Smith 1935f3b1f45cSBarry Smith Logically Collective on Mat 1936f3b1f45cSBarry Smith 1937f3b1f45cSBarry Smith Input Parameters: 1938f3b1f45cSBarry Smith + mat - the shell matrix 1939f3b1f45cSBarry Smith . f - the function 1940f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 1941f3b1f45cSBarry Smith - ctx - an optional context for the function 1942f3b1f45cSBarry Smith 1943f3b1f45cSBarry Smith Output Parameter: 1944f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 1945f3b1f45cSBarry Smith 1946f3b1f45cSBarry Smith Options Database: 1947f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 1948f3b1f45cSBarry Smith 1949f3b1f45cSBarry Smith Level: advanced 1950f3b1f45cSBarry Smith 195195452b02SPatrick Sanan Fortran Notes: 195295452b02SPatrick Sanan Not supported from Fortran 1953f3b1f45cSBarry Smith 1954db781477SPatrick Sanan .seealso: `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellTestMult()` 1955f3b1f45cSBarry Smith @*/ 1956f3b1f45cSBarry Smith PetscErrorCode MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 1957f3b1f45cSBarry Smith { 1958f3b1f45cSBarry Smith Vec x,y,z; 1959f3b1f45cSBarry Smith PetscInt m,n,M,N; 1960f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 1961f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 196274e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 1963f3b1f45cSBarry Smith 1964f3b1f45cSBarry Smith PetscFunctionBegin; 1965f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 19669566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v)); 19679566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat,&x,&y)); 19689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y,&z)); 19699566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat,&m,&n)); 19709566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat,&M,&N)); 19719566063dSJacob Faibussowitsch PetscCall(MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf)); 19729566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(mf,f,ctx)); 19739566063dSJacob Faibussowitsch PetscCall(MatMFFDSetBase(mf,base,NULL)); 19749566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(mf,MATAIJ,&Dmf)); 19759566063dSJacob Faibussowitsch PetscCall(MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf)); 19769566063dSJacob Faibussowitsch PetscCall(MatComputeOperatorTranspose(mat,MATAIJ,&Dmat)); 1977f3b1f45cSBarry Smith 19789566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff)); 19799566063dSJacob Faibussowitsch PetscCall(MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN)); 19809566063dSJacob Faibussowitsch PetscCall(MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm)); 19819566063dSJacob Faibussowitsch PetscCall(MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm)); 1982f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 1983f3b1f45cSBarry Smith flag = PETSC_FALSE; 1984f3b1f45cSBarry Smith if (v) { 19859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm))); 19869566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view")); 19879566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view")); 19889566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view")); 1989f3b1f45cSBarry Smith } 1990f3b1f45cSBarry Smith } else if (v) { 19919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n")); 1992f3b1f45cSBarry Smith } 1993f3b1f45cSBarry Smith if (flg) *flg = flag; 19949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mf)); 19959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Dmat)); 19969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ddiff)); 19979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Dmf)); 19989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 19999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 20009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 2001f3b1f45cSBarry Smith PetscFunctionReturn(0); 2002f3b1f45cSBarry Smith } 2003f3b1f45cSBarry Smith 2004f3b1f45cSBarry Smith /*@C 20050c0fd78eSBarry Smith MatShellSetOperation - Allows user to set a matrix operation for a shell matrix. 2006e51e0e81SBarry Smith 20073f9fe445SBarry Smith Logically Collective on Mat 2008fee21e36SBarry Smith 2009c7fcc2eaSBarry Smith Input Parameters: 2010c7fcc2eaSBarry Smith + mat - the shell matrix 2011c7fcc2eaSBarry Smith . op - the name of the operation 2012789d8953SBarry Smith - g - the function that provides the operation. 2013c7fcc2eaSBarry Smith 201415091d37SBarry Smith Level: advanced 201515091d37SBarry Smith 2016fae171e0SBarry Smith Usage: 2017e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 2018b77ba244SStefano Zampini $ MatCreateShell(comm,m,n,M,N,ctx,&A); 2019b77ba244SStefano Zampini $ MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 20200b627109SLois Curfman McInnes 2021a62d957aSLois Curfman McInnes Notes: 2022e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 20231c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 2024a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 20251c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 2026a62d957aSLois Curfman McInnes 2027a4d85fd7SBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 2028deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 2029deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 2030deebb3c3SLois Curfman McInnes routines, e.g., 2031a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 2032a62d957aSLois Curfman McInnes 20334aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 20344aa34b0aSBarry Smith nonzero on failure. 20354aa34b0aSBarry Smith 2036a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 2037a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 2038a62d957aSLois Curfman McInnes set by MatCreateShell(). 2039a62d957aSLois Curfman McInnes 2040b77ba244SStefano Zampini Use MatSetOperation() to set an operation for any matrix type. For matrix product operations (i.e. MatMat, MatTransposeMat etc) use MatShellSetMatProductOperation() 2041b77ba244SStefano Zampini 204295452b02SPatrick Sanan Fortran Notes: 204395452b02SPatrick Sanan For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 2044c4762a1bSJed Brown generate a matrix. See src/mat/tests/ex120f.F 2045501d9185SBarry Smith 2046db781477SPatrick Sanan .seealso: `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`, `MatShellSetContext()`, `MatSetOperation()`, `MatShellSetManageScalingShifts()`, `MatShellSetMatProductOperation()` 2047e51e0e81SBarry Smith @*/ 2048789d8953SBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*g)(void)) 2049e51e0e81SBarry Smith { 20503a40ed3dSBarry Smith PetscFunctionBegin; 20510700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2052cac4c232SBarry Smith PetscTryMethod(mat,"MatShellSetOperation_C",(Mat,MatOperation,void (*)(void)),(mat,op,g)); 20533a40ed3dSBarry Smith PetscFunctionReturn(0); 2054e51e0e81SBarry Smith } 2055f0479e8cSBarry Smith 2056d4bb536fSBarry Smith /*@C 2057d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 2058d4bb536fSBarry Smith 2059c7fcc2eaSBarry Smith Not Collective 2060c7fcc2eaSBarry Smith 2061d4bb536fSBarry Smith Input Parameters: 2062c7fcc2eaSBarry Smith + mat - the shell matrix 2063c7fcc2eaSBarry Smith - op - the name of the operation 2064d4bb536fSBarry Smith 2065d4bb536fSBarry Smith Output Parameter: 2066789d8953SBarry Smith . g - the function that provides the operation. 2067d4bb536fSBarry Smith 206815091d37SBarry Smith Level: advanced 206915091d37SBarry Smith 2070d4bb536fSBarry Smith Notes: 2071e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 2072d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 2073d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 2074d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 2075d4bb536fSBarry Smith 2076d4bb536fSBarry Smith All user-provided functions have the same calling 2077d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 2078d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 2079d4bb536fSBarry Smith routines, e.g., 2080d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 2081d4bb536fSBarry Smith 2082d4bb536fSBarry Smith Within each user-defined routine, the user should call 2083d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 2084d4bb536fSBarry Smith set by MatCreateShell(). 2085d4bb536fSBarry Smith 2086db781477SPatrick Sanan .seealso: `MatCreateShell()`, `MatShellGetContext()`, `MatShellSetOperation()`, `MatShellSetContext()` 2087d4bb536fSBarry Smith @*/ 2088789d8953SBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**g)(void)) 2089d4bb536fSBarry Smith { 20903a40ed3dSBarry Smith PetscFunctionBegin; 20910700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2092cac4c232SBarry Smith PetscUseMethod(mat,"MatShellGetOperation_C",(Mat,MatOperation,void (**)(void)),(mat,op,g)); 20933a40ed3dSBarry Smith PetscFunctionReturn(0); 2094d4bb536fSBarry Smith } 2095b77ba244SStefano Zampini 2096b77ba244SStefano Zampini /*@ 2097b77ba244SStefano Zampini MatIsShell - Inquires if a matrix is derived from MATSHELL 2098b77ba244SStefano Zampini 2099b77ba244SStefano Zampini Input Parameter: 2100b77ba244SStefano Zampini . mat - the matrix 2101b77ba244SStefano Zampini 2102b77ba244SStefano Zampini Output Parameter: 2103b77ba244SStefano Zampini . flg - the boolean value 2104b77ba244SStefano Zampini 2105b77ba244SStefano Zampini Level: developer 2106b77ba244SStefano Zampini 2107b77ba244SStefano 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) 2108b77ba244SStefano Zampini 2109db781477SPatrick Sanan .seealso: `MatCreateShell()` 2110b77ba244SStefano Zampini @*/ 2111b77ba244SStefano Zampini PetscErrorCode MatIsShell(Mat mat, PetscBool *flg) 2112b77ba244SStefano Zampini { 2113b77ba244SStefano Zampini PetscFunctionBegin; 2114b77ba244SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2115dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg,2); 2116b77ba244SStefano Zampini *flg = (PetscBool)(mat->ops->destroy == MatDestroy_Shell); 2117b77ba244SStefano Zampini PetscFunctionReturn(0); 2118b77ba244SStefano Zampini } 2119