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 */ 239f137db4SBarry Smith PetscErrorCode (*axpy)(Mat,PetscScalar,Mat,MatStructure); 2428f357e3SAlex Fikl PetscErrorCode (*copy)(Mat,Mat,MatStructure); 2528f357e3SAlex Fikl /* 44 */ 266464bf51SAlex Fikl PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode); 2728f357e3SAlex Fikl /* 49 */ 2828f357e3SAlex Fikl /* 54 */ 2928f357e3SAlex Fikl /* 59 */ 3028f357e3SAlex Fikl PetscErrorCode (*destroy)(Mat); 3128f357e3SAlex Fikl /* 64 */ 3228f357e3SAlex Fikl /* 69 */ 3328f357e3SAlex Fikl /* 74 */ 3428f357e3SAlex Fikl /* 79 */ 3528f357e3SAlex Fikl /* 84 */ 3628f357e3SAlex Fikl /* 89 */ 3728f357e3SAlex Fikl /* 94 */ 3828f357e3SAlex Fikl /* 99 */ 3928f357e3SAlex Fikl /* 104 */ 4028f357e3SAlex Fikl /* 109 */ 4128f357e3SAlex Fikl /* 114 */ 4228f357e3SAlex Fikl /* 119 */ 4328f357e3SAlex Fikl /* 124 */ 4428f357e3SAlex Fikl /* 129 */ 4528f357e3SAlex Fikl /* 134 */ 4628f357e3SAlex Fikl /* 139 */ 4728f357e3SAlex Fikl /* 144 */ 4828f357e3SAlex Fikl }; 4928f357e3SAlex Fikl 5028f357e3SAlex Fikl typedef struct { 5128f357e3SAlex Fikl struct _MatShellOps ops[1]; 522205254eSKarl Rupp 53ef66eb69SBarry Smith PetscScalar vscale,vshift; 548fe8eb27SJed Brown Vec dshift; 558fe8eb27SJed Brown Vec left,right; 568fe8eb27SJed Brown Vec left_work,right_work; 575edf6516SJed Brown Vec left_add_work,right_add_work; 589f137db4SBarry Smith Mat axpy; 599f137db4SBarry Smith PetscScalar axpy_vscale; 600c0fd78eSBarry Smith PetscBool managescalingshifts; /* The user will manage the scaling and shifts for the MATSHELL, not the default */ 6120563c6bSBarry Smith void *ctx; 6288cf3e7dSBarry Smith } Mat_Shell; 630c0fd78eSBarry Smith 648fe8eb27SJed Brown /* 650c0fd78eSBarry Smith xx = diag(left)*x 668fe8eb27SJed Brown */ 678fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 688fe8eb27SJed Brown { 698fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 708fe8eb27SJed Brown PetscErrorCode ierr; 718fe8eb27SJed Brown 728fe8eb27SJed Brown PetscFunctionBegin; 730298fd71SBarry Smith *xx = NULL; 748fe8eb27SJed Brown if (!shell->left) { 758fe8eb27SJed Brown *xx = x; 768fe8eb27SJed Brown } else { 778fe8eb27SJed Brown if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);} 788fe8eb27SJed Brown ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr); 798fe8eb27SJed Brown *xx = shell->left_work; 808fe8eb27SJed Brown } 818fe8eb27SJed Brown PetscFunctionReturn(0); 828fe8eb27SJed Brown } 838fe8eb27SJed Brown 840c0fd78eSBarry Smith /* 850c0fd78eSBarry Smith xx = diag(right)*x 860c0fd78eSBarry Smith */ 878fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 888fe8eb27SJed Brown { 898fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 908fe8eb27SJed Brown PetscErrorCode ierr; 918fe8eb27SJed Brown 928fe8eb27SJed Brown PetscFunctionBegin; 930298fd71SBarry Smith *xx = NULL; 948fe8eb27SJed Brown if (!shell->right) { 958fe8eb27SJed Brown *xx = x; 968fe8eb27SJed Brown } else { 978fe8eb27SJed Brown if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);} 988fe8eb27SJed Brown ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr); 998fe8eb27SJed Brown *xx = shell->right_work; 1008fe8eb27SJed Brown } 1018fe8eb27SJed Brown PetscFunctionReturn(0); 1028fe8eb27SJed Brown } 1038fe8eb27SJed Brown 1040c0fd78eSBarry Smith /* 1050c0fd78eSBarry Smith x = diag(left)*x 1060c0fd78eSBarry Smith */ 1078fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 1088fe8eb27SJed Brown { 1098fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1108fe8eb27SJed Brown PetscErrorCode ierr; 1118fe8eb27SJed Brown 1128fe8eb27SJed Brown PetscFunctionBegin; 1138fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);} 1148fe8eb27SJed Brown PetscFunctionReturn(0); 1158fe8eb27SJed Brown } 1168fe8eb27SJed Brown 1170c0fd78eSBarry Smith /* 1180c0fd78eSBarry Smith x = diag(right)*x 1190c0fd78eSBarry Smith */ 1208fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 1218fe8eb27SJed Brown { 1228fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1238fe8eb27SJed Brown PetscErrorCode ierr; 1248fe8eb27SJed Brown 1258fe8eb27SJed Brown PetscFunctionBegin; 1268fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);} 1278fe8eb27SJed Brown PetscFunctionReturn(0); 1288fe8eb27SJed Brown } 1298fe8eb27SJed Brown 1300c0fd78eSBarry Smith /* 1310c0fd78eSBarry Smith Y = vscale*Y + diag(dshift)*X + vshift*X 1320c0fd78eSBarry Smith 1330c0fd78eSBarry Smith On input Y already contains A*x 1340c0fd78eSBarry Smith */ 1358fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 1368fe8eb27SJed Brown { 1378fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1388fe8eb27SJed Brown PetscErrorCode ierr; 1398fe8eb27SJed Brown 1408fe8eb27SJed Brown PetscFunctionBegin; 1418fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 1428fe8eb27SJed Brown PetscInt i,m; 1438fe8eb27SJed Brown const PetscScalar *x,*d; 1448fe8eb27SJed Brown PetscScalar *y; 1458fe8eb27SJed Brown ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr); 1468fe8eb27SJed Brown ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1478fe8eb27SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 1488fe8eb27SJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1498fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 1508fe8eb27SJed Brown ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1518fe8eb27SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 1528fe8eb27SJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 1530c0fd78eSBarry Smith } else { 154027c5fb5SJed Brown ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr); 1558fe8eb27SJed Brown } 156d4c7de66SBarry Smith if (shell->vshift != 0.0) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */ 1578fe8eb27SJed Brown PetscFunctionReturn(0); 1588fe8eb27SJed Brown } 159e51e0e81SBarry Smith 1609d225801SJed Brown /*@ 161a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 162b4fd4287SBarry Smith 163c7fcc2eaSBarry Smith Not Collective 164c7fcc2eaSBarry Smith 165b4fd4287SBarry Smith Input Parameter: 166b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 167b4fd4287SBarry Smith 168b4fd4287SBarry Smith Output Parameter: 169b4fd4287SBarry Smith . ctx - the user provided context 170b4fd4287SBarry Smith 17115091d37SBarry Smith Level: advanced 17215091d37SBarry Smith 173daf670e6SBarry Smith Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this 174daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 175a62d957aSLois Curfman McInnes 176a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context 177a62d957aSLois Curfman McInnes 178ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 1790bc0a719SBarry Smith @*/ 1808fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 181b4fd4287SBarry Smith { 182dfbe8321SBarry Smith PetscErrorCode ierr; 183ace3abfcSBarry Smith PetscBool flg; 184273d9f13SBarry Smith 1853a40ed3dSBarry Smith PetscFunctionBegin; 1860700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1874482741eSBarry Smith PetscValidPointer(ctx,2); 188251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 189940183f3SBarry Smith if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 190ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix"); 1913a40ed3dSBarry Smith PetscFunctionReturn(0); 192b4fd4287SBarry Smith } 193b4fd4287SBarry Smith 194dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 195e51e0e81SBarry Smith { 196dfbe8321SBarry Smith PetscErrorCode ierr; 197bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 198ed3cc1f0SBarry Smith 1993a40ed3dSBarry Smith PetscFunctionBegin; 20028f357e3SAlex Fikl if (shell->ops->destroy) { 20128f357e3SAlex Fikl ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr); 202bf0cc555SLisandro Dalcin } 2030c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 2040c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 2050c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 2068fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 2078fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 2085edf6516SJed Brown ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr); 2095edf6516SJed Brown ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr); 2109f137db4SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 211bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 2123a40ed3dSBarry Smith PetscFunctionReturn(0); 213e51e0e81SBarry Smith } 214e51e0e81SBarry Smith 2157fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str) 2167fabda3fSAlex Fikl { 21728f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data; 2187fabda3fSAlex Fikl PetscErrorCode ierr; 219cb8c8a77SBarry Smith PetscBool matflg; 2207fabda3fSAlex Fikl 2217fabda3fSAlex Fikl PetscFunctionBegin; 22228f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr); 223cb8c8a77SBarry Smith if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL"); 2247fabda3fSAlex Fikl 225cb8c8a77SBarry Smith ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 226cb8c8a77SBarry Smith ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr); 2277fabda3fSAlex Fikl 228cb8c8a77SBarry Smith if (shellA->ops->copy) { 22928f357e3SAlex Fikl ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr); 230cb8c8a77SBarry Smith } 2317fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 2327fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 2330c0fd78eSBarry Smith if (shellA->dshift) { 2340c0fd78eSBarry Smith if (!shellB->dshift) { 2350c0fd78eSBarry Smith ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr); 2367fabda3fSAlex Fikl } 2370c0fd78eSBarry Smith ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr); 2387fabda3fSAlex Fikl } else { 2399f137db4SBarry Smith ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr); 2407fabda3fSAlex Fikl } 2410c0fd78eSBarry Smith if (shellA->left) { 2420c0fd78eSBarry Smith if (!shellB->left) { 2430c0fd78eSBarry Smith ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr); 2447fabda3fSAlex Fikl } 2450c0fd78eSBarry Smith ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr); 2467fabda3fSAlex Fikl } else { 2479f137db4SBarry Smith ierr = VecDestroy(&shellB->left);CHKERRQ(ierr); 2487fabda3fSAlex Fikl } 2490c0fd78eSBarry Smith if (shellA->right) { 2500c0fd78eSBarry Smith if (!shellB->right) { 2510c0fd78eSBarry Smith ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr); 2527fabda3fSAlex Fikl } 2530c0fd78eSBarry Smith ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr); 2547fabda3fSAlex Fikl } else { 2559f137db4SBarry Smith ierr = VecDestroy(&shellB->right);CHKERRQ(ierr); 2567fabda3fSAlex Fikl } 257*40e381c3SBarry Smith ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr); 258*40e381c3SBarry Smith if (shellA->axpy) { 259*40e381c3SBarry Smith ierr = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr); 260*40e381c3SBarry Smith shellB->axpy = shellA->axpy; 261*40e381c3SBarry Smith shellB->axpy_vscale = shellA->axpy_vscale; 262*40e381c3SBarry Smith } 2637fabda3fSAlex Fikl PetscFunctionReturn(0); 2647fabda3fSAlex Fikl } 2657fabda3fSAlex Fikl 266cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M) 267cb8c8a77SBarry Smith { 268cb8c8a77SBarry Smith PetscErrorCode ierr; 269cb8c8a77SBarry Smith void *ctx; 270cb8c8a77SBarry Smith 271cb8c8a77SBarry Smith PetscFunctionBegin; 272cb8c8a77SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 273cb8c8a77SBarry Smith ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr); 274cb8c8a77SBarry Smith ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 275cb8c8a77SBarry Smith PetscFunctionReturn(0); 276cb8c8a77SBarry Smith } 277cb8c8a77SBarry Smith 278dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 279ef66eb69SBarry Smith { 280ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 281dfbe8321SBarry Smith PetscErrorCode ierr; 28225578ef6SJed Brown Vec xx; 283e3f487b0SBarry Smith PetscObjectState instate,outstate; 284ef66eb69SBarry Smith 285ef66eb69SBarry Smith PetscFunctionBegin; 2868fe8eb27SJed Brown ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr); 287e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 28828f357e3SAlex Fikl ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr); 289e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 290e3f487b0SBarry Smith if (instate == outstate) { 291e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 292e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 293e3f487b0SBarry Smith } 2948fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2958fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 2969f137db4SBarry Smith 2979f137db4SBarry Smith if (shell->axpy) { 2989f137db4SBarry Smith if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);} 2999f137db4SBarry Smith ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr); 3009f137db4SBarry Smith ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 3019f137db4SBarry Smith } 302ef66eb69SBarry Smith PetscFunctionReturn(0); 303ef66eb69SBarry Smith } 304ef66eb69SBarry Smith 3055edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 3065edf6516SJed Brown { 3075edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 3085edf6516SJed Brown PetscErrorCode ierr; 3095edf6516SJed Brown 3105edf6516SJed Brown PetscFunctionBegin; 3115edf6516SJed Brown if (y == z) { 3125edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 3135edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 314b8430d49SBarry Smith ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr); 3155edf6516SJed Brown } else { 3165edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 3175edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 3185edf6516SJed Brown } 3195edf6516SJed Brown PetscFunctionReturn(0); 3205edf6516SJed Brown } 3215edf6516SJed Brown 32218be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 32318be62a5SSatish Balay { 32418be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 32518be62a5SSatish Balay PetscErrorCode ierr; 32625578ef6SJed Brown Vec xx; 327e3f487b0SBarry Smith PetscObjectState instate,outstate; 32818be62a5SSatish Balay 32918be62a5SSatish Balay PetscFunctionBegin; 3308fe8eb27SJed Brown ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr); 331e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 3320c0fd78eSBarry Smith if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL"); 33328f357e3SAlex Fikl ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr); 334e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 335e3f487b0SBarry Smith if (instate == outstate) { 336e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 337e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 338e3f487b0SBarry Smith } 3398fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 3408fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 34118be62a5SSatish Balay PetscFunctionReturn(0); 34218be62a5SSatish Balay } 34318be62a5SSatish Balay 3445edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 3455edf6516SJed Brown { 3465edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 3475edf6516SJed Brown PetscErrorCode ierr; 3485edf6516SJed Brown 3495edf6516SJed Brown PetscFunctionBegin; 3505edf6516SJed Brown if (y == z) { 3515edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 3525edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 3535edf6516SJed Brown ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr); 3545edf6516SJed Brown } else { 3555edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 3565edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 3575edf6516SJed Brown } 3585edf6516SJed Brown PetscFunctionReturn(0); 3595edf6516SJed Brown } 3605edf6516SJed Brown 3610c0fd78eSBarry Smith /* 3620c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 3630c0fd78eSBarry Smith */ 36418be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 36518be62a5SSatish Balay { 36618be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 36718be62a5SSatish Balay PetscErrorCode ierr; 36818be62a5SSatish Balay 36918be62a5SSatish Balay PetscFunctionBegin; 3700c0fd78eSBarry Smith if (shell->ops->getdiagonal) { 37128f357e3SAlex Fikl ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr); 3720c0fd78eSBarry Smith } else { 3730c0fd78eSBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 3740c0fd78eSBarry Smith } 37518be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 3768fe8eb27SJed Brown if (shell->dshift) { 3770c0fd78eSBarry Smith ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr); 37818be62a5SSatish Balay } 3790c0fd78eSBarry Smith ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 3808fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 3818fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 38281c02519SBarry Smith if (shell->axpy) { 38381c02519SBarry Smith if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);} 38481c02519SBarry Smith ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr); 38581c02519SBarry Smith ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 38681c02519SBarry Smith } 38718be62a5SSatish Balay PetscFunctionReturn(0); 38818be62a5SSatish Balay } 38918be62a5SSatish Balay 390f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 391ef66eb69SBarry Smith { 392ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3938fe8eb27SJed Brown PetscErrorCode ierr; 394b24ad042SBarry Smith 395ef66eb69SBarry Smith PetscFunctionBegin; 3960c0fd78eSBarry Smith if (shell->left || shell->right) { 3978fe8eb27SJed Brown if (!shell->dshift) { 3980c0fd78eSBarry Smith ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr); 3990c0fd78eSBarry Smith ierr = VecSet(shell->dshift,a);CHKERRQ(ierr); 4000c0fd78eSBarry Smith } else { 4010c0fd78eSBarry Smith if (shell->left) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 4020c0fd78eSBarry Smith if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 4039f137db4SBarry Smith ierr = VecShift(shell->dshift,a);CHKERRQ(ierr); 4040c0fd78eSBarry Smith } 4058fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 4068fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 4078fe8eb27SJed Brown } else shell->vshift += a; 408ef66eb69SBarry Smith PetscFunctionReturn(0); 409ef66eb69SBarry Smith } 410ef66eb69SBarry Smith 4116464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 4126464bf51SAlex Fikl { 4136464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 4146464bf51SAlex Fikl PetscErrorCode ierr; 4156464bf51SAlex Fikl 4166464bf51SAlex Fikl PetscFunctionBegin; 4170c0fd78eSBarry Smith if (ins == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported with INSERT_VALUES"); 4180c0fd78eSBarry Smith if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);} 4190c0fd78eSBarry Smith if (shell->left || shell->right) { 42092fabd1fSBarry Smith if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);} 4219f137db4SBarry Smith if (shell->left && shell->right) { 4229f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 4239f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr); 4249f137db4SBarry Smith } else if (shell->left) { 4259f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 4269f137db4SBarry Smith } else { 4279f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr); 4289f137db4SBarry Smith } 4299f137db4SBarry Smith if (!shell->dshift) { 4309f137db4SBarry Smith ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr); 4319f137db4SBarry Smith ierr = VecCopy(shell->dshift,shell->right_work);CHKERRQ(ierr); 4329f137db4SBarry Smith } else { 4339f137db4SBarry Smith ierr = VecAXPY(shell->dshift,1.0,shell->right_work);CHKERRQ(ierr); 4349f137db4SBarry Smith } 4350c0fd78eSBarry Smith } else { 4360c0fd78eSBarry Smith ierr = VecAXPY(shell->dshift,1.0,D);CHKERRQ(ierr); 4376464bf51SAlex Fikl } 4386464bf51SAlex Fikl PetscFunctionReturn(0); 4396464bf51SAlex Fikl } 4406464bf51SAlex Fikl 441f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 442ef66eb69SBarry Smith { 443ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 4448fe8eb27SJed Brown PetscErrorCode ierr; 445b24ad042SBarry Smith 446ef66eb69SBarry Smith PetscFunctionBegin; 447f4df32b1SMatthew Knepley shell->vscale *= a; 4480c0fd78eSBarry Smith shell->vshift *= a; 4492205254eSKarl Rupp if (shell->dshift) { 4502205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 4510c0fd78eSBarry Smith } 45281c02519SBarry Smith shell->axpy_vscale *= a; 4538fe8eb27SJed Brown PetscFunctionReturn(0); 45418be62a5SSatish Balay } 4558fe8eb27SJed Brown 4568fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 4578fe8eb27SJed Brown { 4588fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 4598fe8eb27SJed Brown PetscErrorCode ierr; 4608fe8eb27SJed Brown 4618fe8eb27SJed Brown PetscFunctionBegin; 46281c02519SBarry Smith if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation"); 4638fe8eb27SJed Brown if (left) { 4640c0fd78eSBarry Smith if (!shell->left) { 4650c0fd78eSBarry Smith ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr); 4668fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 4670c0fd78eSBarry Smith } else { 4680c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 46918be62a5SSatish Balay } 470ef66eb69SBarry Smith } 4718fe8eb27SJed Brown if (right) { 4720c0fd78eSBarry Smith if (!shell->right) { 4730c0fd78eSBarry Smith ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr); 4748fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 4750c0fd78eSBarry Smith } else { 4760c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 4778fe8eb27SJed Brown } 4788fe8eb27SJed Brown } 479ef66eb69SBarry Smith PetscFunctionReturn(0); 480ef66eb69SBarry Smith } 481ef66eb69SBarry Smith 482dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 483ef66eb69SBarry Smith { 484ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 4850c0fd78eSBarry Smith PetscErrorCode ierr; 486ef66eb69SBarry Smith 487ef66eb69SBarry Smith PetscFunctionBegin; 4888fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 489ef66eb69SBarry Smith shell->vshift = 0.0; 490ef66eb69SBarry Smith shell->vscale = 1.0; 4910c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 4920c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 4930c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 49481c02519SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 495ef66eb69SBarry Smith } 496ef66eb69SBarry Smith PetscFunctionReturn(0); 497ef66eb69SBarry Smith } 498ef66eb69SBarry Smith 499cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*); 500b951964fSBarry Smith 5013b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 5023b49f96aSBarry Smith { 5033b49f96aSBarry Smith PetscFunctionBegin; 5043b49f96aSBarry Smith *missing = PETSC_FALSE; 5053b49f96aSBarry Smith PetscFunctionReturn(0); 5063b49f96aSBarry Smith } 5073b49f96aSBarry Smith 5089f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str) 5099f137db4SBarry Smith { 5109f137db4SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 5119f137db4SBarry Smith PetscErrorCode ierr; 5129f137db4SBarry Smith 5139f137db4SBarry Smith PetscFunctionBegin; 5149f137db4SBarry Smith ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 5159f137db4SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 5169f137db4SBarry Smith shell->axpy = X; 5179f137db4SBarry Smith shell->axpy_vscale = a; 5189f137db4SBarry Smith PetscFunctionReturn(0); 5199f137db4SBarry Smith } 5209f137db4SBarry Smith 52109dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 52220563c6bSBarry Smith 0, 52320563c6bSBarry Smith 0, 5249f137db4SBarry Smith 0, 5250c0fd78eSBarry Smith /* 4*/ MatMultAdd_Shell, 5269f137db4SBarry Smith 0, 5270c0fd78eSBarry Smith MatMultTransposeAdd_Shell, 528b951964fSBarry Smith 0, 529b951964fSBarry Smith 0, 530b951964fSBarry Smith 0, 53197304618SKris Buschelman /*10*/ 0, 532b951964fSBarry Smith 0, 533b951964fSBarry Smith 0, 534b951964fSBarry Smith 0, 535b951964fSBarry Smith 0, 53697304618SKris Buschelman /*15*/ 0, 537b951964fSBarry Smith 0, 5380c0fd78eSBarry Smith MatGetDiagonal_Shell, 5398fe8eb27SJed Brown MatDiagonalScale_Shell, 540b951964fSBarry Smith 0, 54197304618SKris Buschelman /*20*/ 0, 542ef66eb69SBarry Smith MatAssemblyEnd_Shell, 543b951964fSBarry Smith 0, 544b951964fSBarry Smith 0, 545d519adbfSMatthew Knepley /*24*/ 0, 546b951964fSBarry Smith 0, 547b951964fSBarry Smith 0, 548b951964fSBarry Smith 0, 549b951964fSBarry Smith 0, 550d519adbfSMatthew Knepley /*29*/ 0, 551b951964fSBarry Smith 0, 552273d9f13SBarry Smith 0, 553b951964fSBarry Smith 0, 554b951964fSBarry Smith 0, 555cb8c8a77SBarry Smith /*34*/ MatDuplicate_Shell, 556b951964fSBarry Smith 0, 557b951964fSBarry Smith 0, 55809dc0095SBarry Smith 0, 55909dc0095SBarry Smith 0, 5609f137db4SBarry Smith /*39*/ MatAXPY_Shell, 56109dc0095SBarry Smith 0, 56209dc0095SBarry Smith 0, 56309dc0095SBarry Smith 0, 564cb8c8a77SBarry Smith MatCopy_Shell, 565d519adbfSMatthew Knepley /*44*/ 0, 566ef66eb69SBarry Smith MatScale_Shell, 567ef66eb69SBarry Smith MatShift_Shell, 5686464bf51SAlex Fikl MatDiagonalSet_Shell, 56909dc0095SBarry Smith 0, 570f73d5cc4SBarry Smith /*49*/ 0, 57109dc0095SBarry Smith 0, 57209dc0095SBarry Smith 0, 57309dc0095SBarry Smith 0, 57409dc0095SBarry Smith 0, 575d519adbfSMatthew Knepley /*54*/ 0, 57609dc0095SBarry Smith 0, 57709dc0095SBarry Smith 0, 57809dc0095SBarry Smith 0, 57909dc0095SBarry Smith 0, 580d519adbfSMatthew Knepley /*59*/ 0, 581b9b97703SBarry Smith MatDestroy_Shell, 58209dc0095SBarry Smith 0, 583357abbc8SBarry Smith 0, 584273d9f13SBarry Smith 0, 585d519adbfSMatthew Knepley /*64*/ 0, 586273d9f13SBarry Smith 0, 587273d9f13SBarry Smith 0, 588273d9f13SBarry Smith 0, 589273d9f13SBarry Smith 0, 590d519adbfSMatthew Knepley /*69*/ 0, 591273d9f13SBarry Smith 0, 592c87e5d42SMatthew Knepley MatConvert_Shell, 593273d9f13SBarry Smith 0, 59497304618SKris Buschelman 0, 595d519adbfSMatthew Knepley /*74*/ 0, 59697304618SKris Buschelman 0, 59797304618SKris Buschelman 0, 59897304618SKris Buschelman 0, 59997304618SKris Buschelman 0, 600d519adbfSMatthew Knepley /*79*/ 0, 60197304618SKris Buschelman 0, 60297304618SKris Buschelman 0, 60397304618SKris Buschelman 0, 604865e5f61SKris Buschelman 0, 605d519adbfSMatthew Knepley /*84*/ 0, 606865e5f61SKris Buschelman 0, 607865e5f61SKris Buschelman 0, 608865e5f61SKris Buschelman 0, 609865e5f61SKris Buschelman 0, 610d519adbfSMatthew Knepley /*89*/ 0, 611865e5f61SKris Buschelman 0, 612865e5f61SKris Buschelman 0, 613865e5f61SKris Buschelman 0, 614865e5f61SKris Buschelman 0, 615d519adbfSMatthew Knepley /*94*/ 0, 616865e5f61SKris Buschelman 0, 617865e5f61SKris Buschelman 0, 6183964eb88SJed Brown 0, 6193964eb88SJed Brown 0, 6203964eb88SJed Brown /*99*/ 0, 6213964eb88SJed Brown 0, 6223964eb88SJed Brown 0, 6233964eb88SJed Brown 0, 6243964eb88SJed Brown 0, 6253964eb88SJed Brown /*104*/ 0, 6263964eb88SJed Brown 0, 6273964eb88SJed Brown 0, 6283964eb88SJed Brown 0, 6293964eb88SJed Brown 0, 6303964eb88SJed Brown /*109*/ 0, 6313964eb88SJed Brown 0, 6323964eb88SJed Brown 0, 6333964eb88SJed Brown 0, 6343b49f96aSBarry Smith MatMissingDiagonal_Shell, 6353964eb88SJed Brown /*114*/ 0, 6363964eb88SJed Brown 0, 6373964eb88SJed Brown 0, 6383964eb88SJed Brown 0, 6393964eb88SJed Brown 0, 6403964eb88SJed Brown /*119*/ 0, 6413964eb88SJed Brown 0, 6423964eb88SJed Brown 0, 6433964eb88SJed Brown 0, 6443964eb88SJed Brown 0, 6453964eb88SJed Brown /*124*/ 0, 6463964eb88SJed Brown 0, 6473964eb88SJed Brown 0, 6483964eb88SJed Brown 0, 6493964eb88SJed Brown 0, 6503964eb88SJed Brown /*129*/ 0, 6513964eb88SJed Brown 0, 6523964eb88SJed Brown 0, 6533964eb88SJed Brown 0, 6543964eb88SJed Brown 0, 6553964eb88SJed Brown /*134*/ 0, 6563964eb88SJed Brown 0, 6573964eb88SJed Brown 0, 6583964eb88SJed Brown 0, 6593964eb88SJed Brown 0, 6603964eb88SJed Brown /*139*/ 0, 661f9426fe0SMark Adams 0, 6623964eb88SJed Brown 0 6633964eb88SJed Brown }; 664273d9f13SBarry Smith 6650bad9183SKris Buschelman /*MC 666fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 6670bad9183SKris Buschelman 6680bad9183SKris Buschelman Level: advanced 6690bad9183SKris Buschelman 6700c0fd78eSBarry Smith .seealso: MatCreateShell() 6710bad9183SKris Buschelman M*/ 6720bad9183SKris Buschelman 6738cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 674273d9f13SBarry Smith { 675273d9f13SBarry Smith Mat_Shell *b; 676dfbe8321SBarry Smith PetscErrorCode ierr; 677273d9f13SBarry Smith 678273d9f13SBarry Smith PetscFunctionBegin; 679273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 680273d9f13SBarry Smith 681b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 682273d9f13SBarry Smith A->data = (void*)b; 683273d9f13SBarry Smith 68426283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 68526283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 686273d9f13SBarry Smith 687273d9f13SBarry Smith b->ctx = 0; 688ef66eb69SBarry Smith b->vshift = 0.0; 689ef66eb69SBarry Smith b->vscale = 1.0; 6900c0fd78eSBarry Smith b->managescalingshifts = PETSC_TRUE; 691273d9f13SBarry Smith A->assembled = PETSC_TRUE; 692210f0121SBarry Smith A->preallocated = PETSC_FALSE; 6932205254eSKarl Rupp 69417667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 695273d9f13SBarry Smith PetscFunctionReturn(0); 696273d9f13SBarry Smith } 697e51e0e81SBarry Smith 6984b828684SBarry Smith /*@C 699052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 700ff756334SLois Curfman McInnes private data storage format. 701e51e0e81SBarry Smith 702c7fcc2eaSBarry Smith Collective on MPI_Comm 703c7fcc2eaSBarry Smith 704e51e0e81SBarry Smith Input Parameters: 705c7fcc2eaSBarry Smith + comm - MPI communicator 706c7fcc2eaSBarry Smith . m - number of local rows (must be given) 707c7fcc2eaSBarry Smith . n - number of local columns (must be given) 708c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 709c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 710c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 711e51e0e81SBarry Smith 712ff756334SLois Curfman McInnes Output Parameter: 71344cd7ae7SLois Curfman McInnes . A - the matrix 714e51e0e81SBarry Smith 715ff2fd236SBarry Smith Level: advanced 716ff2fd236SBarry Smith 717f39d1f56SLois Curfman McInnes Usage: 7187b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 719f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 720c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 721f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 722f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 723f39d1f56SLois Curfman McInnes 724ff756334SLois Curfman McInnes Notes: 725ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 726ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 727ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 728e51e0e81SBarry Smith 729daf670e6SBarry Smith Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this 730daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 731daf670e6SBarry Smith in as the ctx argument. 732f38a8d56SBarry Smith 733f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 734f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 735645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 736645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 737645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 738645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 739645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 740645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 741645985a0SLois Curfman McInnes For example, 742f39d1f56SLois Curfman McInnes 7430c0fd78eSBarry Smith MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), and MatScale() internally so these 7440c0fd78eSBarry Smith operations cannot be overwritten unless MatShellSetManageScalingShifts() is called. 7450c0fd78eSBarry Smith 746f39d1f56SLois Curfman McInnes $ 747f39d1f56SLois Curfman McInnes $ Vec x, y 7487b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 749645985a0SLois Curfman McInnes $ Mat A 750f39d1f56SLois Curfman McInnes $ 751c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 752c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 753f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 754c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 755c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 756c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 757645985a0SLois Curfman McInnes $ MatMult(A,x,y); 758645985a0SLois Curfman McInnes $ MatDestroy(A); 759f39d1f56SLois Curfman McInnes $ VecDestroy(y); VecDestroy(x); 760645985a0SLois Curfman McInnes $ 761e51e0e81SBarry Smith 7620c0fd78eSBarry Smith For rectangular matrices do all the scalings and shifts make sense? 7630c0fd78eSBarry Smith 7640c0fd78eSBarry Smith Developers Notes: Regarding shifting and scaling. The general form is 7650c0fd78eSBarry Smith 7660c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 7670c0fd78eSBarry Smith 7680c0fd78eSBarry Smith The order you apply the operations is important. For example if you have a dshift then 7690c0fd78eSBarry Smith apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift 7700c0fd78eSBarry Smith you get s*vscale*A + diag(shift) 7710c0fd78eSBarry Smith 7720c0fd78eSBarry Smith A is the user provided function. 7730c0fd78eSBarry Smith 7740b627109SLois Curfman McInnes .keywords: matrix, shell, create 7750b627109SLois Curfman McInnes 7760c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts() 777e51e0e81SBarry Smith @*/ 7787087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 779e51e0e81SBarry Smith { 780dfbe8321SBarry Smith PetscErrorCode ierr; 781ed3cc1f0SBarry Smith 7823a40ed3dSBarry Smith PetscFunctionBegin; 783f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 784f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 785273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 786273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 7870fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 788273d9f13SBarry Smith PetscFunctionReturn(0); 789c7fcc2eaSBarry Smith } 790c7fcc2eaSBarry Smith 791c6866cfdSSatish Balay /*@ 792273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 793c7fcc2eaSBarry Smith 7943f9fe445SBarry Smith Logically Collective on Mat 795c7fcc2eaSBarry Smith 796273d9f13SBarry Smith Input Parameters: 797273d9f13SBarry Smith + mat - the shell matrix 798273d9f13SBarry Smith - ctx - the context 799273d9f13SBarry Smith 800273d9f13SBarry Smith Level: advanced 801273d9f13SBarry Smith 802daf670e6SBarry Smith Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this 803daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 804273d9f13SBarry Smith 805273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 8060bc0a719SBarry Smith @*/ 8077087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 808273d9f13SBarry Smith { 809273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 810dfbe8321SBarry Smith PetscErrorCode ierr; 811ace3abfcSBarry Smith PetscBool flg; 812273d9f13SBarry Smith 813273d9f13SBarry Smith PetscFunctionBegin; 8140700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 815251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 816273d9f13SBarry Smith if (flg) { 817273d9f13SBarry Smith shell->ctx = ctx; 818ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix"); 8193a40ed3dSBarry Smith PetscFunctionReturn(0); 820e51e0e81SBarry Smith } 821e51e0e81SBarry Smith 8220c0fd78eSBarry Smith /*@ 8230c0fd78eSBarry Smith MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately 8240c0fd78eSBarry Smith after MatCreateShell() 8250c0fd78eSBarry Smith 8260c0fd78eSBarry Smith Logically Collective on Mat 8270c0fd78eSBarry Smith 8280c0fd78eSBarry Smith Input Parameter: 8290c0fd78eSBarry Smith . mat - the shell matrix 8300c0fd78eSBarry Smith 8310c0fd78eSBarry Smith Level: advanced 8320c0fd78eSBarry Smith 8330c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation() 8340c0fd78eSBarry Smith @*/ 8350c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A) 8360c0fd78eSBarry Smith { 8370c0fd78eSBarry Smith PetscErrorCode ierr; 8380c0fd78eSBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 8390c0fd78eSBarry Smith PetscBool flg; 8400c0fd78eSBarry Smith 8410c0fd78eSBarry Smith PetscFunctionBegin; 8420c0fd78eSBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 8430c0fd78eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr); 8440c0fd78eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 8450c0fd78eSBarry Smith shell->managescalingshifts = PETSC_FALSE; 8460c0fd78eSBarry Smith A->ops->diagonalscale = 0; 8470c0fd78eSBarry Smith A->ops->scale = 0; 8480c0fd78eSBarry Smith A->ops->shift = 0; 8490c0fd78eSBarry Smith A->ops->diagonalset = 0; 8500c0fd78eSBarry Smith PetscFunctionReturn(0); 8510c0fd78eSBarry Smith } 8520c0fd78eSBarry Smith 853c16cb8f2SBarry Smith /*@C 854f3b1f45cSBarry Smith MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function. 855f3b1f45cSBarry Smith 856f3b1f45cSBarry Smith Logically Collective on Mat 857f3b1f45cSBarry Smith 858f3b1f45cSBarry Smith Input Parameters: 859f3b1f45cSBarry Smith + mat - the shell matrix 860f3b1f45cSBarry Smith . f - the function 861f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 862f3b1f45cSBarry Smith - ctx - an optional context for the function 863f3b1f45cSBarry Smith 864f3b1f45cSBarry Smith Output Parameter: 865f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 866f3b1f45cSBarry Smith 867f3b1f45cSBarry Smith Options Database: 868f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 869f3b1f45cSBarry Smith 870f3b1f45cSBarry Smith Level: advanced 871f3b1f45cSBarry Smith 872f3b1f45cSBarry Smith Fortran Notes: Not supported from Fortran 873f3b1f45cSBarry Smith 874f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose() 875f3b1f45cSBarry Smith @*/ 876f3b1f45cSBarry Smith PetscErrorCode MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 877f3b1f45cSBarry Smith { 878f3b1f45cSBarry Smith PetscErrorCode ierr; 879f3b1f45cSBarry Smith PetscInt m,n; 880f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 881f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 88274e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 883f3b1f45cSBarry Smith 884f3b1f45cSBarry Smith PetscFunctionBegin; 885f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 886f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr); 887f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 888f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr); 889f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 890f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 891f3b1f45cSBarry Smith 892f3b1f45cSBarry Smith ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr); 893f3b1f45cSBarry Smith ierr = MatComputeExplicitOperator(mat,&Dmat);CHKERRQ(ierr); 894f3b1f45cSBarry Smith 895f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 896f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 897f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 898f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 899f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 900f3b1f45cSBarry Smith flag = PETSC_FALSE; 901f3b1f45cSBarry Smith if (v) { 902f3b1f45cSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm)); 903f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 904f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 905f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 906f3b1f45cSBarry Smith } 907f3b1f45cSBarry Smith } else if (v) { 908f3b1f45cSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n"); 909f3b1f45cSBarry Smith } 910f3b1f45cSBarry Smith if (flg) *flg = flag; 911f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 912f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 913f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 914f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 915f3b1f45cSBarry Smith PetscFunctionReturn(0); 916f3b1f45cSBarry Smith } 917f3b1f45cSBarry Smith 918f3b1f45cSBarry Smith /*@C 919f3b1f45cSBarry Smith MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function. 920f3b1f45cSBarry Smith 921f3b1f45cSBarry Smith Logically Collective on Mat 922f3b1f45cSBarry Smith 923f3b1f45cSBarry Smith Input Parameters: 924f3b1f45cSBarry Smith + mat - the shell matrix 925f3b1f45cSBarry Smith . f - the function 926f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 927f3b1f45cSBarry Smith - ctx - an optional context for the function 928f3b1f45cSBarry Smith 929f3b1f45cSBarry Smith Output Parameter: 930f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 931f3b1f45cSBarry Smith 932f3b1f45cSBarry Smith Options Database: 933f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 934f3b1f45cSBarry Smith 935f3b1f45cSBarry Smith Level: advanced 936f3b1f45cSBarry Smith 937f3b1f45cSBarry Smith Fortran Notes: Not supported from Fortran 938f3b1f45cSBarry Smith 939f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult() 940f3b1f45cSBarry Smith @*/ 941f3b1f45cSBarry Smith PetscErrorCode MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 942f3b1f45cSBarry Smith { 943f3b1f45cSBarry Smith PetscErrorCode ierr; 944f3b1f45cSBarry Smith Vec x,y,z; 945f3b1f45cSBarry Smith PetscInt m,n,M,N; 946f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 947f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 94874e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 949f3b1f45cSBarry Smith 950f3b1f45cSBarry Smith PetscFunctionBegin; 951f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 952f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr); 953f3b1f45cSBarry Smith ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr); 954f3b1f45cSBarry Smith ierr = VecDuplicate(y,&z);CHKERRQ(ierr); 955f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 956f3b1f45cSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 957f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr); 958f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 959f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 960f3b1f45cSBarry Smith ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr); 961f3b1f45cSBarry Smith ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr); 962f3b1f45cSBarry Smith ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr); 963f3b1f45cSBarry Smith 964f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 965f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 966f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 967f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 968f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 969f3b1f45cSBarry Smith flag = PETSC_FALSE; 970f3b1f45cSBarry Smith if (v) { 971f3b1f45cSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm)); 972f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 973f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 974f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 975f3b1f45cSBarry Smith } 976f3b1f45cSBarry Smith } else if (v) { 977f3b1f45cSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n"); 978f3b1f45cSBarry Smith } 979f3b1f45cSBarry Smith if (flg) *flg = flag; 980f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 981f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 982f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 983f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 984f3b1f45cSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 985f3b1f45cSBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 986f3b1f45cSBarry Smith ierr = VecDestroy(&z);CHKERRQ(ierr); 987f3b1f45cSBarry Smith PetscFunctionReturn(0); 988f3b1f45cSBarry Smith } 989f3b1f45cSBarry Smith 990f3b1f45cSBarry Smith /*@C 9913a3eedf2SBarry Smith MatShellSetOperation - Allows user to set a matrix operation for 9923a3eedf2SBarry Smith a shell matrix. 9930c0fd78eSBarry Smith MatShellSetOperation - Allows user to set a matrix operation for a shell matrix. 994e51e0e81SBarry Smith 9953f9fe445SBarry Smith Logically Collective on Mat 996fee21e36SBarry Smith 997c7fcc2eaSBarry Smith Input Parameters: 998c7fcc2eaSBarry Smith + mat - the shell matrix 999c7fcc2eaSBarry Smith . op - the name of the operation 1000c7fcc2eaSBarry Smith - f - the function that provides the operation. 1001c7fcc2eaSBarry Smith 100215091d37SBarry Smith Level: advanced 100315091d37SBarry Smith 1004fae171e0SBarry Smith Usage: 1005e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 1006f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 1007c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 10080b627109SLois Curfman McInnes 1009a62d957aSLois Curfman McInnes Notes: 1010e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 10111c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 1012a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 10131c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 1014a62d957aSLois Curfman McInnes 1015a4d85fd7SBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 1016deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 1017deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 1018deebb3c3SLois Curfman McInnes routines, e.g., 1019a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1020a62d957aSLois Curfman McInnes 10214aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 10224aa34b0aSBarry Smith nonzero on failure. 10234aa34b0aSBarry Smith 1024a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 1025a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 1026a62d957aSLois Curfman McInnes set by MatCreateShell(). 1027a62d957aSLois Curfman McInnes 10282a7a6963SBarry Smith Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 1029501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 1030501d9185SBarry Smith 10310c0fd78eSBarry Smith Use MatSetOperation() to set an operation for any matrix type 10320c0fd78eSBarry Smith 1033a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation 1034a62d957aSLois Curfman McInnes 10350c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts() 1036e51e0e81SBarry Smith @*/ 10377087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 1038e51e0e81SBarry Smith { 1039dfbe8321SBarry Smith PetscErrorCode ierr; 1040ace3abfcSBarry Smith PetscBool flg; 10410c0fd78eSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 1042273d9f13SBarry Smith 10433a40ed3dSBarry Smith PetscFunctionBegin; 10440700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 10450c0fd78eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 10460c0fd78eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 10475edf6516SJed Brown switch (op) { 10485edf6516SJed Brown case MATOP_DESTROY: 104928f357e3SAlex Fikl shell->ops->destroy = (PetscErrorCode (*)(Mat))f; 10505edf6516SJed Brown break; 10516464bf51SAlex Fikl case MATOP_DIAGONAL_SET: 10520c0fd78eSBarry Smith case MATOP_DIAGONAL_SCALE: 10530c0fd78eSBarry Smith case MATOP_SHIFT: 10540c0fd78eSBarry Smith case MATOP_SCALE: 10559f137db4SBarry Smith case MATOP_AXPY: 10560c0fd78eSBarry Smith if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()"); 10570c0fd78eSBarry Smith (((void(**)(void))mat->ops)[op]) = f; 10586464bf51SAlex Fikl break; 10590c0fd78eSBarry Smith case MATOP_GET_DIAGONAL: 10600c0fd78eSBarry Smith if (shell->managescalingshifts) shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 10610c0fd78eSBarry Smith else mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 10625edf6516SJed Brown case MATOP_VIEW: 106340a5a6c4SBarry Smith if (!mat->ops->viewnative) { 106440a5a6c4SBarry Smith mat->ops->viewnative = mat->ops->view; 106540a5a6c4SBarry Smith } 10665edf6516SJed Brown mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 10675edf6516SJed Brown break; 10685edf6516SJed Brown case MATOP_MULT: 10699f137db4SBarry Smith if (shell->managescalingshifts){ 10709f137db4SBarry Smith shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 10719f137db4SBarry Smith mat->ops->mult = MatMult_Shell; 10729f137db4SBarry Smith } else mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 10735edf6516SJed Brown break; 10745edf6516SJed Brown case MATOP_MULT_TRANSPOSE: 10759f137db4SBarry Smith if (shell->managescalingshifts) { 10769f137db4SBarry Smith shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 10775259c5a4SBarry Smith mat->ops->multtranspose = MatMultTranspose_Shell; 10789f137db4SBarry Smith } else mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 10795edf6516SJed Brown break; 10805edf6516SJed Brown default: 10815edf6516SJed Brown (((void(**)(void))mat->ops)[op]) = f; 1082a62d957aSLois Curfman McInnes } 10833a40ed3dSBarry Smith PetscFunctionReturn(0); 1084e51e0e81SBarry Smith } 1085f0479e8cSBarry Smith 1086d4bb536fSBarry Smith /*@C 1087d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 1088d4bb536fSBarry Smith 1089c7fcc2eaSBarry Smith Not Collective 1090c7fcc2eaSBarry Smith 1091d4bb536fSBarry Smith Input Parameters: 1092c7fcc2eaSBarry Smith + mat - the shell matrix 1093c7fcc2eaSBarry Smith - op - the name of the operation 1094d4bb536fSBarry Smith 1095d4bb536fSBarry Smith Output Parameter: 1096d4bb536fSBarry Smith . f - the function that provides the operation. 1097d4bb536fSBarry Smith 109815091d37SBarry Smith Level: advanced 109915091d37SBarry Smith 1100d4bb536fSBarry Smith Notes: 1101e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 1102d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 1103d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 1104d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 1105d4bb536fSBarry Smith 1106d4bb536fSBarry Smith All user-provided functions have the same calling 1107d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 1108d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 1109d4bb536fSBarry Smith routines, e.g., 1110d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1111d4bb536fSBarry Smith 1112d4bb536fSBarry Smith Within each user-defined routine, the user should call 1113d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 1114d4bb536fSBarry Smith set by MatCreateShell(). 1115d4bb536fSBarry Smith 1116d4bb536fSBarry Smith .keywords: matrix, shell, set, operation 1117d4bb536fSBarry Smith 1118ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 1119d4bb536fSBarry Smith @*/ 11207087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 1121d4bb536fSBarry Smith { 1122dfbe8321SBarry Smith PetscErrorCode ierr; 1123ace3abfcSBarry Smith PetscBool flg; 1124273d9f13SBarry Smith 11253a40ed3dSBarry Smith PetscFunctionBegin; 11260700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 112728f357e3SAlex Fikl switch (op) { 112828f357e3SAlex Fikl case MATOP_DESTROY: 1129251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 1130273d9f13SBarry Smith if (flg) { 1131d4bb536fSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 113228f357e3SAlex Fikl *f = (void (*)(void))shell->ops->destroy; 1133c7fcc2eaSBarry Smith } else { 1134c134de8dSSatish Balay *f = (void (*)(void))mat->ops->destroy; 1135d4bb536fSBarry Smith } 113628f357e3SAlex Fikl break; 113728f357e3SAlex Fikl case MATOP_DIAGONAL_SET: 113828f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 113928f357e3SAlex Fikl if (flg) { 114028f357e3SAlex Fikl Mat_Shell *shell = (Mat_Shell*)mat->data; 114128f357e3SAlex Fikl *f = (void (*)(void))shell->ops->diagonalset; 1142c7fcc2eaSBarry Smith } else { 114328f357e3SAlex Fikl *f = (void (*)(void))mat->ops->diagonalset; 114428f357e3SAlex Fikl } 114528f357e3SAlex Fikl break; 114628f357e3SAlex Fikl case MATOP_VIEW: 114728f357e3SAlex Fikl *f = (void (*)(void))mat->ops->view; 114828f357e3SAlex Fikl break; 114928f357e3SAlex Fikl default: 1150c134de8dSSatish Balay *f = (((void (**)(void))mat->ops)[op]); 1151d4bb536fSBarry Smith } 11523a40ed3dSBarry Smith PetscFunctionReturn(0); 1153d4bb536fSBarry Smith } 1154d4bb536fSBarry Smith 11550c0fd78eSBarry Smith /*@C 11560c0fd78eSBarry Smith MatSetOperation - Allows user to set a matrix operation for any matrix type 11570c0fd78eSBarry Smith 11580c0fd78eSBarry Smith Logically Collective on Mat 11590c0fd78eSBarry Smith 11600c0fd78eSBarry Smith Input Parameters: 11610c0fd78eSBarry Smith + mat - the shell matrix 11620c0fd78eSBarry Smith . op - the name of the operation 11630c0fd78eSBarry Smith - f - the function that provides the operation. 11640c0fd78eSBarry Smith 11650c0fd78eSBarry Smith Level: developer 11660c0fd78eSBarry Smith 11670c0fd78eSBarry Smith Usage: 11680c0fd78eSBarry Smith $ extern PetscErrorCode usermult(Mat,Vec,Vec); 11690c0fd78eSBarry Smith $ ierr = MatCreateXXX(comm,...&A); 11700c0fd78eSBarry Smith $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 11710c0fd78eSBarry Smith 11720c0fd78eSBarry Smith Notes: 11730c0fd78eSBarry Smith See the file include/petscmat.h for a complete list of matrix 11740c0fd78eSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 11750c0fd78eSBarry Smith <OPERATION> is the name (in all capital letters) of the 11760c0fd78eSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 11770c0fd78eSBarry Smith 11780c0fd78eSBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 11790c0fd78eSBarry Smith sequence as the usual matrix interface routines, since they 11800c0fd78eSBarry Smith are intended to be accessed via the usual matrix interface 11810c0fd78eSBarry Smith routines, e.g., 11820c0fd78eSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 11830c0fd78eSBarry Smith 11840c0fd78eSBarry Smith In particular each function MUST return an error code of 0 on success and 11850c0fd78eSBarry Smith nonzero on failure. 11860c0fd78eSBarry Smith 11870c0fd78eSBarry Smith This routine is distinct from MatShellSetOperation() in that it can be called on any matrix type. 11880c0fd78eSBarry Smith 11890c0fd78eSBarry Smith .keywords: matrix, shell, set, operation 11900c0fd78eSBarry Smith 11910c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation() 11920c0fd78eSBarry Smith @*/ 11930c0fd78eSBarry Smith PetscErrorCode MatSetOperation(Mat mat,MatOperation op,void (*f)(void)) 11940c0fd78eSBarry Smith { 11950c0fd78eSBarry Smith PetscFunctionBegin; 11960c0fd78eSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 11970c0fd78eSBarry Smith (((void(**)(void))mat->ops)[op]) = f; 11980c0fd78eSBarry Smith PetscFunctionReturn(0); 11990c0fd78eSBarry Smith } 1200