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*/ 9e51e0e81SBarry Smith 1028f357e3SAlex Fikl struct _MatShellOps { 1128f357e3SAlex Fikl /* 0 */ 126849ba73SBarry Smith PetscErrorCode (*mult)(Mat,Vec,Vec); 1328f357e3SAlex Fikl /* 5 */ 1418be62a5SSatish Balay PetscErrorCode (*multtranspose)(Mat,Vec,Vec); 1528f357e3SAlex Fikl /* 10 */ 1628f357e3SAlex Fikl /* 15 */ 1718be62a5SSatish Balay PetscErrorCode (*getdiagonal)(Mat,Vec); 1828f357e3SAlex Fikl /* 20 */ 1928f357e3SAlex Fikl /* 24 */ 2028f357e3SAlex Fikl /* 29 */ 2128f357e3SAlex Fikl /* 34 */ 2228f357e3SAlex Fikl /* 39 */ 2328f357e3SAlex Fikl PetscErrorCode (*copy)(Mat,Mat,MatStructure); 2428f357e3SAlex Fikl /* 44 */ 256464bf51SAlex Fikl PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode); 2628f357e3SAlex Fikl /* 49 */ 2728f357e3SAlex Fikl /* 54 */ 2828f357e3SAlex Fikl /* 59 */ 2928f357e3SAlex Fikl PetscErrorCode (*destroy)(Mat); 3028f357e3SAlex Fikl /* 64 */ 3128f357e3SAlex Fikl /* 69 */ 3228f357e3SAlex Fikl /* 74 */ 3328f357e3SAlex Fikl /* 79 */ 3428f357e3SAlex Fikl /* 84 */ 3528f357e3SAlex Fikl /* 89 */ 3628f357e3SAlex Fikl /* 94 */ 3728f357e3SAlex Fikl /* 99 */ 3828f357e3SAlex Fikl /* 104 */ 3928f357e3SAlex Fikl /* 109 */ 4028f357e3SAlex Fikl /* 114 */ 4128f357e3SAlex Fikl /* 119 */ 4228f357e3SAlex Fikl /* 124 */ 4328f357e3SAlex Fikl /* 129 */ 4428f357e3SAlex Fikl /* 134 */ 4528f357e3SAlex Fikl /* 139 */ 4628f357e3SAlex Fikl /* 144 */ 4728f357e3SAlex Fikl }; 4828f357e3SAlex Fikl 4928f357e3SAlex Fikl typedef struct { 5028f357e3SAlex Fikl struct _MatShellOps ops[1]; 512205254eSKarl Rupp 52ef66eb69SBarry Smith PetscScalar vscale,vshift; 538fe8eb27SJed Brown Vec dshift; 548fe8eb27SJed Brown Vec left,right; 558fe8eb27SJed Brown Vec dshift_owned,left_owned,right_owned; 568fe8eb27SJed Brown Vec left_work,right_work; 575edf6516SJed Brown Vec left_add_work,right_add_work; 588fe8eb27SJed Brown PetscBool usingscaled; 5920563c6bSBarry Smith void *ctx; 6088cf3e7dSBarry Smith } Mat_Shell; 618fe8eb27SJed Brown /* 628fe8eb27SJed Brown The most general expression for the matrix is 638fe8eb27SJed Brown 648fe8eb27SJed Brown S = L (a A + B) R 658fe8eb27SJed Brown 668fe8eb27SJed Brown where 678fe8eb27SJed Brown A is the matrix defined by the user's function 688fe8eb27SJed Brown a is a scalar multiple 698fe8eb27SJed Brown L is left scaling 708fe8eb27SJed Brown R is right scaling 718fe8eb27SJed Brown B is a diagonal shift defined by 728fe8eb27SJed Brown diag(dshift) if the vector dshift is non-NULL 738fe8eb27SJed Brown vscale*identity otherwise 748fe8eb27SJed Brown 758fe8eb27SJed Brown The following identities apply: 768fe8eb27SJed Brown 778fe8eb27SJed Brown Scale by c: 788fe8eb27SJed Brown c [L (a A + B) R] = L [(a c) A + c B] R 798fe8eb27SJed Brown 808fe8eb27SJed Brown Shift by c: 818fe8eb27SJed Brown [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R 828fe8eb27SJed Brown 838fe8eb27SJed Brown Diagonal scaling is achieved by simply multiplying with existing L and R vectors 848fe8eb27SJed Brown 858fe8eb27SJed Brown In the data structure: 868fe8eb27SJed Brown 878fe8eb27SJed Brown vscale=1.0 means no special scaling will be applied 888fe8eb27SJed Brown dshift=NULL means a constant diagonal shift (fall back to vshift) 898fe8eb27SJed Brown vshift=0.0 means no constant diagonal shift, note that vshift is only used if dshift is NULL 908fe8eb27SJed Brown */ 918fe8eb27SJed Brown 928fe8eb27SJed Brown static PetscErrorCode MatMult_Shell(Mat,Vec,Vec); 938fe8eb27SJed Brown static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec); 948fe8eb27SJed Brown static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec); 9528f357e3SAlex Fikl static PetscErrorCode MatCopy_Shell(Mat,Mat,MatStructure); 968fe8eb27SJed Brown 978fe8eb27SJed Brown static PetscErrorCode MatShellUseScaledMethods(Mat Y) 988fe8eb27SJed Brown { 998fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 1008fe8eb27SJed Brown 1018fe8eb27SJed Brown PetscFunctionBegin; 1028fe8eb27SJed Brown if (shell->usingscaled) PetscFunctionReturn(0); 10328f357e3SAlex Fikl shell->ops->mult = Y->ops->mult; 1048fe8eb27SJed Brown Y->ops->mult = MatMult_Shell; 1058fe8eb27SJed Brown if (Y->ops->multtranspose) { 10628f357e3SAlex Fikl shell->ops->multtranspose = Y->ops->multtranspose; 1078fe8eb27SJed Brown Y->ops->multtranspose = MatMultTranspose_Shell; 1088fe8eb27SJed Brown } 1098fe8eb27SJed Brown if (Y->ops->getdiagonal) { 11028f357e3SAlex Fikl shell->ops->getdiagonal = Y->ops->getdiagonal; 1118fe8eb27SJed Brown Y->ops->getdiagonal = MatGetDiagonal_Shell; 1128fe8eb27SJed Brown } 11328f357e3SAlex Fikl if (Y->ops->copy) { 11428f357e3SAlex Fikl shell->ops->copy = Y->ops->copy; 11528f357e3SAlex Fikl Y->ops->copy = MatCopy_Shell; 11628f357e3SAlex Fikl } 1178fe8eb27SJed Brown shell->usingscaled = PETSC_TRUE; 1188fe8eb27SJed Brown PetscFunctionReturn(0); 1198fe8eb27SJed Brown } 1208fe8eb27SJed Brown 1218fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 1228fe8eb27SJed Brown { 1238fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1248fe8eb27SJed Brown PetscErrorCode ierr; 1258fe8eb27SJed Brown 1268fe8eb27SJed Brown PetscFunctionBegin; 1270298fd71SBarry Smith *xx = NULL; 1288fe8eb27SJed Brown if (!shell->left) { 1298fe8eb27SJed Brown *xx = x; 1308fe8eb27SJed Brown } else { 1318fe8eb27SJed Brown if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);} 1328fe8eb27SJed Brown ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr); 1338fe8eb27SJed Brown *xx = shell->left_work; 1348fe8eb27SJed Brown } 1358fe8eb27SJed Brown PetscFunctionReturn(0); 1368fe8eb27SJed Brown } 1378fe8eb27SJed Brown 1388fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 1398fe8eb27SJed Brown { 1408fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1418fe8eb27SJed Brown PetscErrorCode ierr; 1428fe8eb27SJed Brown 1438fe8eb27SJed Brown PetscFunctionBegin; 1440298fd71SBarry Smith *xx = NULL; 1458fe8eb27SJed Brown if (!shell->right) { 1468fe8eb27SJed Brown *xx = x; 1478fe8eb27SJed Brown } else { 1488fe8eb27SJed Brown if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);} 1498fe8eb27SJed Brown ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr); 1508fe8eb27SJed Brown *xx = shell->right_work; 1518fe8eb27SJed Brown } 1528fe8eb27SJed Brown PetscFunctionReturn(0); 1538fe8eb27SJed Brown } 1548fe8eb27SJed Brown 1558fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 1568fe8eb27SJed Brown { 1578fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1588fe8eb27SJed Brown PetscErrorCode ierr; 1598fe8eb27SJed Brown 1608fe8eb27SJed Brown PetscFunctionBegin; 1618fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);} 1628fe8eb27SJed Brown PetscFunctionReturn(0); 1638fe8eb27SJed Brown } 1648fe8eb27SJed Brown 1658fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 1668fe8eb27SJed Brown { 1678fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1688fe8eb27SJed Brown PetscErrorCode ierr; 1698fe8eb27SJed Brown 1708fe8eb27SJed Brown PetscFunctionBegin; 1718fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);} 1728fe8eb27SJed Brown PetscFunctionReturn(0); 1738fe8eb27SJed Brown } 1748fe8eb27SJed Brown 1758fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 1768fe8eb27SJed Brown { 1778fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1788fe8eb27SJed Brown PetscErrorCode ierr; 1798fe8eb27SJed Brown 1808fe8eb27SJed Brown PetscFunctionBegin; 1818fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 1828fe8eb27SJed Brown PetscInt i,m; 1838fe8eb27SJed Brown const PetscScalar *x,*d; 1848fe8eb27SJed Brown PetscScalar *y; 1858fe8eb27SJed Brown ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr); 1868fe8eb27SJed Brown ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1878fe8eb27SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 1888fe8eb27SJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1898fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 1908fe8eb27SJed Brown ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1918fe8eb27SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 1928fe8eb27SJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 193880538faSHong Zhang } else if (PetscAbsScalar(shell->vshift) != 0) { 1948fe8eb27SJed Brown ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr); 195027c5fb5SJed Brown } else if (shell->vscale != 1.0) { 196027c5fb5SJed Brown ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr); 1978fe8eb27SJed Brown } 1988fe8eb27SJed Brown PetscFunctionReturn(0); 1998fe8eb27SJed Brown } 200e51e0e81SBarry Smith 2019d225801SJed Brown /*@ 202a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 203b4fd4287SBarry Smith 204c7fcc2eaSBarry Smith Not Collective 205c7fcc2eaSBarry Smith 206b4fd4287SBarry Smith Input Parameter: 207b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 208b4fd4287SBarry Smith 209b4fd4287SBarry Smith Output Parameter: 210b4fd4287SBarry Smith . ctx - the user provided context 211b4fd4287SBarry Smith 21215091d37SBarry Smith Level: advanced 21315091d37SBarry Smith 214daf670e6SBarry Smith Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this 215daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 216a62d957aSLois Curfman McInnes 217a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context 218a62d957aSLois Curfman McInnes 219ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 2200bc0a719SBarry Smith @*/ 2218fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 222b4fd4287SBarry Smith { 223dfbe8321SBarry Smith PetscErrorCode ierr; 224ace3abfcSBarry Smith PetscBool flg; 225273d9f13SBarry Smith 2263a40ed3dSBarry Smith PetscFunctionBegin; 2270700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2284482741eSBarry Smith PetscValidPointer(ctx,2); 229251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 230940183f3SBarry Smith if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 231ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix"); 2323a40ed3dSBarry Smith PetscFunctionReturn(0); 233b4fd4287SBarry Smith } 234b4fd4287SBarry Smith 235dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 236e51e0e81SBarry Smith { 237dfbe8321SBarry Smith PetscErrorCode ierr; 238bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 239ed3cc1f0SBarry Smith 2403a40ed3dSBarry Smith PetscFunctionBegin; 24128f357e3SAlex Fikl if (shell->ops->destroy) { 24228f357e3SAlex Fikl ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr); 243bf0cc555SLisandro Dalcin } 2448fe8eb27SJed Brown ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr); 2458fe8eb27SJed Brown ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr); 2468fe8eb27SJed Brown ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr); 2478fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 2488fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 2495edf6516SJed Brown ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr); 2505edf6516SJed Brown ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr); 251bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 2523a40ed3dSBarry Smith PetscFunctionReturn(0); 253e51e0e81SBarry Smith } 254e51e0e81SBarry Smith 2557fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str) 2567fabda3fSAlex Fikl { 25728f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data; 25828f357e3SAlex Fikl PetscBool matflg,shellflg; 2597fabda3fSAlex Fikl PetscErrorCode ierr; 2607fabda3fSAlex Fikl 2617fabda3fSAlex Fikl PetscFunctionBegin; 26228f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr); 26328f357e3SAlex Fikl if(!matflg) { SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL"); } 2647fabda3fSAlex Fikl 2657fabda3fSAlex Fikl ierr = MatShellUseScaledMethods(B);CHKERRQ(ierr); 26628f357e3SAlex Fikl ierr = PetscMemcmp(A->ops,B->ops,sizeof(struct _MatOps),&matflg);CHKERRQ(ierr); 26728f357e3SAlex Fikl ierr = PetscMemcmp(shellA->ops,shellB->ops,sizeof(struct _MatShellOps),&shellflg);CHKERRQ(ierr); 26828f357e3SAlex Fikl if (!matflg || !shellflg) { SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Matrices cannot be copied because they have different operations"); } 2697fabda3fSAlex Fikl 27028f357e3SAlex Fikl ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr); 2717fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 2727fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 2737fabda3fSAlex Fikl if (shellA->dshift_owned) { 2747fabda3fSAlex Fikl if (!shellB->dshift_owned) { 2757fabda3fSAlex Fikl ierr = VecDuplicate(shellA->dshift_owned,&shellB->dshift_owned);CHKERRQ(ierr); 2767fabda3fSAlex Fikl } 2777fabda3fSAlex Fikl ierr = VecCopy(shellA->dshift_owned,shellB->dshift_owned);CHKERRQ(ierr); 2787fabda3fSAlex Fikl shellB->dshift = shellB->dshift_owned; 2797fabda3fSAlex Fikl } else { 2807fabda3fSAlex Fikl shellB->dshift = NULL; 2817fabda3fSAlex Fikl } 2827fabda3fSAlex Fikl if (shellA->left_owned) { 2837fabda3fSAlex Fikl if (!shellB->left_owned) { 2847fabda3fSAlex Fikl ierr = VecDuplicate(shellA->left_owned,&shellB->left_owned);CHKERRQ(ierr); 2857fabda3fSAlex Fikl } 2867fabda3fSAlex Fikl ierr = VecCopy(shellA->left_owned,shellB->left_owned);CHKERRQ(ierr); 2877fabda3fSAlex Fikl shellB->left = shellB->left_owned; 2887fabda3fSAlex Fikl } else { 2897fabda3fSAlex Fikl shellB->left = NULL; 2907fabda3fSAlex Fikl } 2917fabda3fSAlex Fikl if (shellA->right_owned) { 2927fabda3fSAlex Fikl if (!shellB->right_owned) { 2937fabda3fSAlex Fikl ierr = VecDuplicate(shellA->right_owned,&shellB->right_owned);CHKERRQ(ierr); 2947fabda3fSAlex Fikl } 2957fabda3fSAlex Fikl ierr = VecCopy(shellA->right_owned,shellB->right_owned);CHKERRQ(ierr); 2967fabda3fSAlex Fikl shellB->right = shellB->right_owned; 2977fabda3fSAlex Fikl } else { 2987fabda3fSAlex Fikl shellB->right = NULL; 2997fabda3fSAlex Fikl } 3007fabda3fSAlex Fikl PetscFunctionReturn(0); 3017fabda3fSAlex Fikl } 3027fabda3fSAlex Fikl 303dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 304ef66eb69SBarry Smith { 305ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 306dfbe8321SBarry Smith PetscErrorCode ierr; 30725578ef6SJed Brown Vec xx; 308e3f487b0SBarry Smith PetscObjectState instate,outstate; 309ef66eb69SBarry Smith 310ef66eb69SBarry Smith PetscFunctionBegin; 3118fe8eb27SJed Brown ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr); 312e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 31328f357e3SAlex Fikl ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr); 314e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 315e3f487b0SBarry Smith if (instate == outstate) { 316e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 317e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 318e3f487b0SBarry Smith } 3198fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 3208fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 321ef66eb69SBarry Smith PetscFunctionReturn(0); 322ef66eb69SBarry Smith } 323ef66eb69SBarry Smith 3245edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 3255edf6516SJed Brown { 3265edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 3275edf6516SJed Brown PetscErrorCode ierr; 3285edf6516SJed Brown 3295edf6516SJed Brown PetscFunctionBegin; 3305edf6516SJed Brown if (y == z) { 3315edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 3325edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 333b8430d49SBarry Smith ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr); 3345edf6516SJed Brown } else { 3355edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 3365edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 3375edf6516SJed Brown } 3385edf6516SJed Brown PetscFunctionReturn(0); 3395edf6516SJed Brown } 3405edf6516SJed Brown 34118be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 34218be62a5SSatish Balay { 34318be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 34418be62a5SSatish Balay PetscErrorCode ierr; 34525578ef6SJed Brown Vec xx; 346e3f487b0SBarry Smith PetscObjectState instate,outstate; 34718be62a5SSatish Balay 34818be62a5SSatish Balay PetscFunctionBegin; 3498fe8eb27SJed Brown ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr); 350e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 35128f357e3SAlex Fikl ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr); 352e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 353e3f487b0SBarry Smith if (instate == outstate) { 354e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 355e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 356e3f487b0SBarry Smith } 3578fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 3588fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 35918be62a5SSatish Balay PetscFunctionReturn(0); 36018be62a5SSatish Balay } 36118be62a5SSatish Balay 3625edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 3635edf6516SJed Brown { 3645edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 3655edf6516SJed Brown PetscErrorCode ierr; 3665edf6516SJed Brown 3675edf6516SJed Brown PetscFunctionBegin; 3685edf6516SJed Brown if (y == z) { 3695edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 3705edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 3715edf6516SJed Brown ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr); 3725edf6516SJed Brown } else { 3735edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 3745edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 3755edf6516SJed Brown } 3765edf6516SJed Brown PetscFunctionReturn(0); 3775edf6516SJed Brown } 3785edf6516SJed Brown 37918be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 38018be62a5SSatish Balay { 38118be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 38218be62a5SSatish Balay PetscErrorCode ierr; 38318be62a5SSatish Balay 38418be62a5SSatish Balay PetscFunctionBegin; 38528f357e3SAlex Fikl ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr); 38618be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 3878fe8eb27SJed Brown if (shell->dshift) { 3888fe8eb27SJed Brown ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr); 3898fe8eb27SJed Brown } else { 39018be62a5SSatish Balay ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 39118be62a5SSatish Balay } 3928fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 3938fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 39418be62a5SSatish Balay PetscFunctionReturn(0); 39518be62a5SSatish Balay } 39618be62a5SSatish Balay 397f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 398ef66eb69SBarry Smith { 399ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 4008fe8eb27SJed Brown PetscErrorCode ierr; 401b24ad042SBarry Smith 402ef66eb69SBarry Smith PetscFunctionBegin; 4038fe8eb27SJed Brown if (shell->left || shell->right || shell->dshift) { 4048fe8eb27SJed Brown if (!shell->dshift) { 4058fe8eb27SJed Brown if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);} 4068fe8eb27SJed Brown shell->dshift = shell->dshift_owned; 4078fe8eb27SJed Brown ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr); 4088fe8eb27SJed Brown } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);} 4098fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 4108fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 4118fe8eb27SJed Brown } else shell->vshift += a; 4128fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 413ef66eb69SBarry Smith PetscFunctionReturn(0); 414ef66eb69SBarry Smith } 415ef66eb69SBarry Smith 4166464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 4176464bf51SAlex Fikl { 4186464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 4196464bf51SAlex Fikl PetscErrorCode ierr; 4206464bf51SAlex Fikl 4216464bf51SAlex Fikl PetscFunctionBegin; 4226464bf51SAlex Fikl if (shell->vscale != 1.0 || shell->left || shell->right) { 4236464bf51SAlex Fikl SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported on scaled matrix"); 4246464bf51SAlex Fikl } 4256464bf51SAlex Fikl 42628f357e3SAlex Fikl if (shell->ops->diagonalset) { 42728f357e3SAlex Fikl ierr = (*shell->ops->diagonalset)(A,D,ins);CHKERRQ(ierr); 42828f357e3SAlex Fikl } 4296464bf51SAlex Fikl shell->vshift = 0.0; 4306464bf51SAlex Fikl if (shell->dshift) { ierr = VecZeroEntries(shell->dshift);CHKERRQ(ierr); } 4316464bf51SAlex Fikl PetscFunctionReturn(0); 4326464bf51SAlex Fikl } 4336464bf51SAlex Fikl 434f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 435ef66eb69SBarry Smith { 436ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 4378fe8eb27SJed Brown PetscErrorCode ierr; 438b24ad042SBarry Smith 439ef66eb69SBarry Smith PetscFunctionBegin; 440f4df32b1SMatthew Knepley shell->vscale *= a; 4412205254eSKarl Rupp if (shell->dshift) { 4422205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 4432205254eSKarl Rupp } else shell->vshift *= a; 4448fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 4458fe8eb27SJed Brown PetscFunctionReturn(0); 44618be62a5SSatish Balay } 4478fe8eb27SJed Brown 4488fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 4498fe8eb27SJed Brown { 4508fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 4518fe8eb27SJed Brown PetscErrorCode ierr; 4528fe8eb27SJed Brown 4538fe8eb27SJed Brown PetscFunctionBegin; 4548fe8eb27SJed Brown if (left) { 4558fe8eb27SJed Brown if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);} 4562205254eSKarl Rupp if (shell->left) { 4572205254eSKarl Rupp ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 4582205254eSKarl Rupp } else { 4598fe8eb27SJed Brown shell->left = shell->left_owned; 4608fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 46118be62a5SSatish Balay } 462ef66eb69SBarry Smith } 4638fe8eb27SJed Brown if (right) { 4648fe8eb27SJed Brown if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);} 4652205254eSKarl Rupp if (shell->right) { 4662205254eSKarl Rupp ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 4672205254eSKarl Rupp } else { 4688fe8eb27SJed Brown shell->right = shell->right_owned; 4698fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 4708fe8eb27SJed Brown } 4718fe8eb27SJed Brown } 4728fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 473ef66eb69SBarry Smith PetscFunctionReturn(0); 474ef66eb69SBarry Smith } 475ef66eb69SBarry Smith 476dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 477ef66eb69SBarry Smith { 478ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 479ef66eb69SBarry Smith 480ef66eb69SBarry Smith PetscFunctionBegin; 4818fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 482ef66eb69SBarry Smith shell->vshift = 0.0; 483ef66eb69SBarry Smith shell->vscale = 1.0; 4840298fd71SBarry Smith shell->dshift = NULL; 4850298fd71SBarry Smith shell->left = NULL; 4860298fd71SBarry Smith shell->right = NULL; 48728f357e3SAlex Fikl if (shell->ops->mult) { 48828f357e3SAlex Fikl Y->ops->mult = shell->ops->mult; 48928f357e3SAlex Fikl shell->ops->mult = NULL; 4907fb7d96aSJed Brown } 49128f357e3SAlex Fikl if (shell->ops->multtranspose) { 49228f357e3SAlex Fikl Y->ops->multtranspose = shell->ops->multtranspose; 49328f357e3SAlex Fikl shell->ops->multtranspose = NULL; 4947fb7d96aSJed Brown } 49528f357e3SAlex Fikl if (shell->ops->getdiagonal) { 49628f357e3SAlex Fikl Y->ops->getdiagonal = shell->ops->getdiagonal; 49728f357e3SAlex Fikl shell->ops->getdiagonal = NULL; 49828f357e3SAlex Fikl } 49928f357e3SAlex Fikl if (shell->ops->copy) { 50028f357e3SAlex Fikl Y->ops->copy = shell->ops->copy; 50128f357e3SAlex Fikl shell->ops->copy = NULL; 5027fb7d96aSJed Brown } 5037fb7d96aSJed Brown shell->usingscaled = PETSC_FALSE; 504ef66eb69SBarry Smith } 505ef66eb69SBarry Smith PetscFunctionReturn(0); 506ef66eb69SBarry Smith } 507ef66eb69SBarry Smith 508cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*); 509b951964fSBarry Smith 5103b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 5113b49f96aSBarry Smith { 5123b49f96aSBarry Smith PetscFunctionBegin; 5133b49f96aSBarry Smith *missing = PETSC_FALSE; 5143b49f96aSBarry Smith PetscFunctionReturn(0); 5153b49f96aSBarry Smith } 5163b49f96aSBarry Smith 51709dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 51820563c6bSBarry Smith 0, 51920563c6bSBarry Smith 0, 52020563c6bSBarry Smith 0, 52197304618SKris Buschelman /* 4*/ 0, 52220563c6bSBarry Smith 0, 523b951964fSBarry Smith 0, 524b951964fSBarry Smith 0, 525b951964fSBarry Smith 0, 526b951964fSBarry Smith 0, 52797304618SKris Buschelman /*10*/ 0, 528b951964fSBarry Smith 0, 529b951964fSBarry Smith 0, 530b951964fSBarry Smith 0, 531b951964fSBarry Smith 0, 53297304618SKris Buschelman /*15*/ 0, 533b951964fSBarry Smith 0, 534b951964fSBarry Smith 0, 5358fe8eb27SJed Brown MatDiagonalScale_Shell, 536b951964fSBarry Smith 0, 53797304618SKris Buschelman /*20*/ 0, 538ef66eb69SBarry Smith MatAssemblyEnd_Shell, 539b951964fSBarry Smith 0, 540b951964fSBarry Smith 0, 541d519adbfSMatthew Knepley /*24*/ 0, 542b951964fSBarry Smith 0, 543b951964fSBarry Smith 0, 544b951964fSBarry Smith 0, 545b951964fSBarry Smith 0, 546d519adbfSMatthew Knepley /*29*/ 0, 547b951964fSBarry Smith 0, 548273d9f13SBarry Smith 0, 549b951964fSBarry Smith 0, 550b951964fSBarry Smith 0, 551d519adbfSMatthew Knepley /*34*/ 0, 552b951964fSBarry Smith 0, 553b951964fSBarry Smith 0, 55409dc0095SBarry Smith 0, 55509dc0095SBarry Smith 0, 556d519adbfSMatthew Knepley /*39*/ 0, 55709dc0095SBarry Smith 0, 55809dc0095SBarry Smith 0, 55909dc0095SBarry Smith 0, 56028f357e3SAlex Fikl 0, 561d519adbfSMatthew Knepley /*44*/ 0, 562ef66eb69SBarry Smith MatScale_Shell, 563ef66eb69SBarry Smith MatShift_Shell, 5646464bf51SAlex Fikl MatDiagonalSet_Shell, 56509dc0095SBarry Smith 0, 566f73d5cc4SBarry Smith /*49*/ 0, 56709dc0095SBarry Smith 0, 56809dc0095SBarry Smith 0, 56909dc0095SBarry Smith 0, 57009dc0095SBarry Smith 0, 571d519adbfSMatthew Knepley /*54*/ 0, 57209dc0095SBarry Smith 0, 57309dc0095SBarry Smith 0, 57409dc0095SBarry Smith 0, 57509dc0095SBarry Smith 0, 576d519adbfSMatthew Knepley /*59*/ 0, 577b9b97703SBarry Smith MatDestroy_Shell, 57809dc0095SBarry Smith 0, 579357abbc8SBarry Smith 0, 580273d9f13SBarry Smith 0, 581d519adbfSMatthew Knepley /*64*/ 0, 582273d9f13SBarry Smith 0, 583273d9f13SBarry Smith 0, 584273d9f13SBarry Smith 0, 585273d9f13SBarry Smith 0, 586d519adbfSMatthew Knepley /*69*/ 0, 587273d9f13SBarry Smith 0, 588c87e5d42SMatthew Knepley MatConvert_Shell, 589273d9f13SBarry Smith 0, 59097304618SKris Buschelman 0, 591d519adbfSMatthew Knepley /*74*/ 0, 59297304618SKris Buschelman 0, 59397304618SKris Buschelman 0, 59497304618SKris Buschelman 0, 59597304618SKris Buschelman 0, 596d519adbfSMatthew Knepley /*79*/ 0, 59797304618SKris Buschelman 0, 59897304618SKris Buschelman 0, 59997304618SKris Buschelman 0, 600865e5f61SKris Buschelman 0, 601d519adbfSMatthew Knepley /*84*/ 0, 602865e5f61SKris Buschelman 0, 603865e5f61SKris Buschelman 0, 604865e5f61SKris Buschelman 0, 605865e5f61SKris Buschelman 0, 606d519adbfSMatthew Knepley /*89*/ 0, 607865e5f61SKris Buschelman 0, 608865e5f61SKris Buschelman 0, 609865e5f61SKris Buschelman 0, 610865e5f61SKris Buschelman 0, 611d519adbfSMatthew Knepley /*94*/ 0, 612865e5f61SKris Buschelman 0, 613865e5f61SKris Buschelman 0, 6143964eb88SJed Brown 0, 6153964eb88SJed Brown 0, 6163964eb88SJed Brown /*99*/ 0, 6173964eb88SJed Brown 0, 6183964eb88SJed Brown 0, 6193964eb88SJed Brown 0, 6203964eb88SJed Brown 0, 6213964eb88SJed Brown /*104*/ 0, 6223964eb88SJed Brown 0, 6233964eb88SJed Brown 0, 6243964eb88SJed Brown 0, 6253964eb88SJed Brown 0, 6263964eb88SJed Brown /*109*/ 0, 6273964eb88SJed Brown 0, 6283964eb88SJed Brown 0, 6293964eb88SJed Brown 0, 6303b49f96aSBarry Smith MatMissingDiagonal_Shell, 6313964eb88SJed Brown /*114*/ 0, 6323964eb88SJed Brown 0, 6333964eb88SJed Brown 0, 6343964eb88SJed Brown 0, 6353964eb88SJed Brown 0, 6363964eb88SJed Brown /*119*/ 0, 6373964eb88SJed Brown 0, 6383964eb88SJed Brown 0, 6393964eb88SJed Brown 0, 6403964eb88SJed Brown 0, 6413964eb88SJed Brown /*124*/ 0, 6423964eb88SJed Brown 0, 6433964eb88SJed Brown 0, 6443964eb88SJed Brown 0, 6453964eb88SJed Brown 0, 6463964eb88SJed Brown /*129*/ 0, 6473964eb88SJed Brown 0, 6483964eb88SJed Brown 0, 6493964eb88SJed Brown 0, 6503964eb88SJed Brown 0, 6513964eb88SJed Brown /*134*/ 0, 6523964eb88SJed Brown 0, 6533964eb88SJed Brown 0, 6543964eb88SJed Brown 0, 6553964eb88SJed Brown 0, 6563964eb88SJed Brown /*139*/ 0, 657f9426fe0SMark Adams 0, 6583964eb88SJed Brown 0 6593964eb88SJed Brown }; 660273d9f13SBarry Smith 6610bad9183SKris Buschelman /*MC 662fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 6630bad9183SKris Buschelman 6640bad9183SKris Buschelman Level: advanced 6650bad9183SKris Buschelman 6660bad9183SKris Buschelman .seealso: MatCreateShell 6670bad9183SKris Buschelman M*/ 6680bad9183SKris Buschelman 6698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 670273d9f13SBarry Smith { 671273d9f13SBarry Smith Mat_Shell *b; 672dfbe8321SBarry Smith PetscErrorCode ierr; 673273d9f13SBarry Smith 674273d9f13SBarry Smith PetscFunctionBegin; 675273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 676273d9f13SBarry Smith 677b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 678273d9f13SBarry Smith A->data = (void*)b; 679273d9f13SBarry Smith 68026283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 68126283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 682273d9f13SBarry Smith 683273d9f13SBarry Smith b->ctx = 0; 684ef66eb69SBarry Smith b->vshift = 0.0; 685ef66eb69SBarry Smith b->vscale = 1.0; 686273d9f13SBarry Smith A->assembled = PETSC_TRUE; 687210f0121SBarry Smith A->preallocated = PETSC_FALSE; 6882205254eSKarl Rupp 68917667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 690273d9f13SBarry Smith PetscFunctionReturn(0); 691273d9f13SBarry Smith } 692e51e0e81SBarry Smith 6934b828684SBarry Smith /*@C 694052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 695ff756334SLois Curfman McInnes private data storage format. 696e51e0e81SBarry Smith 697c7fcc2eaSBarry Smith Collective on MPI_Comm 698c7fcc2eaSBarry Smith 699e51e0e81SBarry Smith Input Parameters: 700c7fcc2eaSBarry Smith + comm - MPI communicator 701c7fcc2eaSBarry Smith . m - number of local rows (must be given) 702c7fcc2eaSBarry Smith . n - number of local columns (must be given) 703c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 704c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 705c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 706e51e0e81SBarry Smith 707ff756334SLois Curfman McInnes Output Parameter: 70844cd7ae7SLois Curfman McInnes . A - the matrix 709e51e0e81SBarry Smith 710ff2fd236SBarry Smith Level: advanced 711ff2fd236SBarry Smith 712f39d1f56SLois Curfman McInnes Usage: 7137b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 714f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 715c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 716f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 717f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 718f39d1f56SLois Curfman McInnes 719ff756334SLois Curfman McInnes Notes: 720ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 721ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 722ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 723e51e0e81SBarry Smith 724daf670e6SBarry Smith Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this 725daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 726daf670e6SBarry Smith in as the ctx argument. 727f38a8d56SBarry Smith 728f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 729f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 730645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 731645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 732645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 733645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 734645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 735645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 736645985a0SLois Curfman McInnes For example, 737f39d1f56SLois Curfman McInnes 738f39d1f56SLois Curfman McInnes $ 739f39d1f56SLois Curfman McInnes $ Vec x, y 7407b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 741645985a0SLois Curfman McInnes $ Mat A 742f39d1f56SLois Curfman McInnes $ 743c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 744c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 745f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 746c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 747c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 748c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 749645985a0SLois Curfman McInnes $ MatMult(A,x,y); 750645985a0SLois Curfman McInnes $ MatDestroy(A); 751f39d1f56SLois Curfman McInnes $ VecDestroy(y); VecDestroy(x); 752645985a0SLois Curfman McInnes $ 753e51e0e81SBarry Smith 7540b627109SLois Curfman McInnes .keywords: matrix, shell, create 7550b627109SLois Curfman McInnes 756ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext() 757e51e0e81SBarry Smith @*/ 7587087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 759e51e0e81SBarry Smith { 760dfbe8321SBarry Smith PetscErrorCode ierr; 761ed3cc1f0SBarry Smith 7623a40ed3dSBarry Smith PetscFunctionBegin; 763f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 764f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 765273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 766273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 7670fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 768273d9f13SBarry Smith PetscFunctionReturn(0); 769c7fcc2eaSBarry Smith } 770c7fcc2eaSBarry Smith 771c6866cfdSSatish Balay /*@ 772273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 773c7fcc2eaSBarry Smith 7743f9fe445SBarry Smith Logically Collective on Mat 775c7fcc2eaSBarry Smith 776273d9f13SBarry Smith Input Parameters: 777273d9f13SBarry Smith + mat - the shell matrix 778273d9f13SBarry Smith - ctx - the context 779273d9f13SBarry Smith 780273d9f13SBarry Smith Level: advanced 781273d9f13SBarry Smith 782daf670e6SBarry Smith Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this 783daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 784273d9f13SBarry Smith 785273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 7860bc0a719SBarry Smith @*/ 7877087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 788273d9f13SBarry Smith { 789273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 790dfbe8321SBarry Smith PetscErrorCode ierr; 791ace3abfcSBarry Smith PetscBool flg; 792273d9f13SBarry Smith 793273d9f13SBarry Smith PetscFunctionBegin; 7940700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 795251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 796273d9f13SBarry Smith if (flg) { 797273d9f13SBarry Smith shell->ctx = ctx; 798ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix"); 7993a40ed3dSBarry Smith PetscFunctionReturn(0); 800e51e0e81SBarry Smith } 801e51e0e81SBarry Smith 802c16cb8f2SBarry Smith /*@C 8033a3eedf2SBarry Smith MatShellSetOperation - Allows user to set a matrix operation for 8043a3eedf2SBarry Smith a shell matrix. 805e51e0e81SBarry Smith 8063f9fe445SBarry Smith Logically Collective on Mat 807fee21e36SBarry Smith 808c7fcc2eaSBarry Smith Input Parameters: 809c7fcc2eaSBarry Smith + mat - the shell matrix 810c7fcc2eaSBarry Smith . op - the name of the operation 811c7fcc2eaSBarry Smith - f - the function that provides the operation. 812c7fcc2eaSBarry Smith 81315091d37SBarry Smith Level: advanced 81415091d37SBarry Smith 815fae171e0SBarry Smith Usage: 816e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 817f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 818c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 8190b627109SLois Curfman McInnes 820a62d957aSLois Curfman McInnes Notes: 821e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 8221c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 823a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 8241c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 825a62d957aSLois Curfman McInnes 82625296bd5SBarry Smith All user-provided functions (execept for MATOP_DESTROY) should have the same calling 827deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 828deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 829deebb3c3SLois Curfman McInnes routines, e.g., 830a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 831a62d957aSLois Curfman McInnes 8324aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 8334aa34b0aSBarry Smith nonzero on failure. 8344aa34b0aSBarry Smith 835a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 836a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 837a62d957aSLois Curfman McInnes set by MatCreateShell(). 838a62d957aSLois Curfman McInnes 8392a7a6963SBarry Smith Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 840501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 841501d9185SBarry Smith 842a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation 843a62d957aSLois Curfman McInnes 844ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext() 845e51e0e81SBarry Smith @*/ 8467087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 847e51e0e81SBarry Smith { 848dfbe8321SBarry Smith PetscErrorCode ierr; 849ace3abfcSBarry Smith PetscBool flg; 850273d9f13SBarry Smith 8513a40ed3dSBarry Smith PetscFunctionBegin; 8520700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 8535edf6516SJed Brown switch (op) { 8545edf6516SJed Brown case MATOP_DESTROY: 855251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 856273d9f13SBarry Smith if (flg) { 857a62d957aSLois Curfman McInnes Mat_Shell *shell = (Mat_Shell*)mat->data; 85828f357e3SAlex Fikl shell->ops->destroy = (PetscErrorCode (*)(Mat))f; 8596849ba73SBarry Smith } else mat->ops->destroy = (PetscErrorCode (*)(Mat))f; 8605edf6516SJed Brown break; 8616464bf51SAlex Fikl case MATOP_DIAGONAL_SET: 8626464bf51SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 8636464bf51SAlex Fikl if (flg) { 8646464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)mat->data; 86528f357e3SAlex Fikl shell->ops->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f; 8666464bf51SAlex Fikl } else { 8676464bf51SAlex Fikl mat->ops->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f; 8686464bf51SAlex Fikl } 8696464bf51SAlex Fikl break; 8705edf6516SJed Brown case MATOP_VIEW: 8715edf6516SJed Brown mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 8725edf6516SJed Brown break; 8735edf6516SJed Brown case MATOP_MULT: 8745edf6516SJed Brown mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 875*875c1aa2SLisandro Dalcin if (!mat->ops->multadd) { 876*875c1aa2SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 877*875c1aa2SLisandro Dalcin if (flg) mat->ops->multadd = MatMultAdd_Shell; 878*875c1aa2SLisandro Dalcin } 8795edf6516SJed Brown break; 8805edf6516SJed Brown case MATOP_MULT_TRANSPOSE: 8815edf6516SJed Brown mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 882*875c1aa2SLisandro Dalcin if (!mat->ops->multtransposeadd) { 883*875c1aa2SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 884*875c1aa2SLisandro Dalcin if (flg) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell; 885*875c1aa2SLisandro Dalcin } 8865edf6516SJed Brown break; 8875edf6516SJed Brown default: 8885edf6516SJed Brown (((void(**)(void))mat->ops)[op]) = f; 889a62d957aSLois Curfman McInnes } 8903a40ed3dSBarry Smith PetscFunctionReturn(0); 891e51e0e81SBarry Smith } 892f0479e8cSBarry Smith 893d4bb536fSBarry Smith /*@C 894d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 895d4bb536fSBarry Smith 896c7fcc2eaSBarry Smith Not Collective 897c7fcc2eaSBarry Smith 898d4bb536fSBarry Smith Input Parameters: 899c7fcc2eaSBarry Smith + mat - the shell matrix 900c7fcc2eaSBarry Smith - op - the name of the operation 901d4bb536fSBarry Smith 902d4bb536fSBarry Smith Output Parameter: 903d4bb536fSBarry Smith . f - the function that provides the operation. 904d4bb536fSBarry Smith 90515091d37SBarry Smith Level: advanced 90615091d37SBarry Smith 907d4bb536fSBarry Smith Notes: 908e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 909d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 910d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 911d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 912d4bb536fSBarry Smith 913d4bb536fSBarry Smith All user-provided functions have the same calling 914d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 915d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 916d4bb536fSBarry Smith routines, e.g., 917d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 918d4bb536fSBarry Smith 919d4bb536fSBarry Smith Within each user-defined routine, the user should call 920d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 921d4bb536fSBarry Smith set by MatCreateShell(). 922d4bb536fSBarry Smith 923d4bb536fSBarry Smith .keywords: matrix, shell, set, operation 924d4bb536fSBarry Smith 925ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 926d4bb536fSBarry Smith @*/ 9277087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 928d4bb536fSBarry Smith { 929dfbe8321SBarry Smith PetscErrorCode ierr; 930ace3abfcSBarry Smith PetscBool flg; 931273d9f13SBarry Smith 9323a40ed3dSBarry Smith PetscFunctionBegin; 9330700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 93428f357e3SAlex Fikl switch (op) { 93528f357e3SAlex Fikl case MATOP_DESTROY: 936251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 937273d9f13SBarry Smith if (flg) { 938d4bb536fSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 93928f357e3SAlex Fikl *f = (void (*)(void))shell->ops->destroy; 940c7fcc2eaSBarry Smith } else { 941c134de8dSSatish Balay *f = (void (*)(void))mat->ops->destroy; 942d4bb536fSBarry Smith } 94328f357e3SAlex Fikl break; 94428f357e3SAlex Fikl case MATOP_DIAGONAL_SET: 94528f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 94628f357e3SAlex Fikl if (flg) { 94728f357e3SAlex Fikl Mat_Shell *shell = (Mat_Shell*)mat->data; 94828f357e3SAlex Fikl *f = (void (*)(void))shell->ops->diagonalset; 949c7fcc2eaSBarry Smith } else { 95028f357e3SAlex Fikl *f = (void (*)(void))mat->ops->diagonalset; 95128f357e3SAlex Fikl } 95228f357e3SAlex Fikl break; 95328f357e3SAlex Fikl case MATOP_VIEW: 95428f357e3SAlex Fikl *f = (void (*)(void))mat->ops->view; 95528f357e3SAlex Fikl break; 95628f357e3SAlex Fikl default: 957c134de8dSSatish Balay *f = (((void (**)(void))mat->ops)[op]); 958d4bb536fSBarry Smith } 9593a40ed3dSBarry Smith PetscFunctionReturn(0); 960d4bb536fSBarry Smith } 961d4bb536fSBarry Smith 962