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 left_work,right_work; 565edf6516SJed Brown Vec left_add_work,right_add_work; 57*0c0fd78eSBarry Smith PetscBool managescalingshifts; /* The user will manage the scaling and shifts for the MATSHELL, not the default */ 5820563c6bSBarry Smith void *ctx; 5988cf3e7dSBarry Smith } Mat_Shell; 60*0c0fd78eSBarry Smith 618fe8eb27SJed Brown /* 62*0c0fd78eSBarry Smith xx = diag(left)*x 638fe8eb27SJed Brown */ 648fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 658fe8eb27SJed Brown { 668fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 678fe8eb27SJed Brown PetscErrorCode ierr; 688fe8eb27SJed Brown 698fe8eb27SJed Brown PetscFunctionBegin; 700298fd71SBarry Smith *xx = NULL; 718fe8eb27SJed Brown if (!shell->left) { 728fe8eb27SJed Brown *xx = x; 738fe8eb27SJed Brown } else { 748fe8eb27SJed Brown if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);} 758fe8eb27SJed Brown ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr); 768fe8eb27SJed Brown *xx = shell->left_work; 778fe8eb27SJed Brown } 788fe8eb27SJed Brown PetscFunctionReturn(0); 798fe8eb27SJed Brown } 808fe8eb27SJed Brown 81*0c0fd78eSBarry Smith /* 82*0c0fd78eSBarry Smith xx = diag(right)*x 83*0c0fd78eSBarry Smith */ 848fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 858fe8eb27SJed Brown { 868fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 878fe8eb27SJed Brown PetscErrorCode ierr; 888fe8eb27SJed Brown 898fe8eb27SJed Brown PetscFunctionBegin; 900298fd71SBarry Smith *xx = NULL; 918fe8eb27SJed Brown if (!shell->right) { 928fe8eb27SJed Brown *xx = x; 938fe8eb27SJed Brown } else { 948fe8eb27SJed Brown if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);} 958fe8eb27SJed Brown ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr); 968fe8eb27SJed Brown *xx = shell->right_work; 978fe8eb27SJed Brown } 988fe8eb27SJed Brown PetscFunctionReturn(0); 998fe8eb27SJed Brown } 1008fe8eb27SJed Brown 101*0c0fd78eSBarry Smith /* 102*0c0fd78eSBarry Smith x = diag(left)*x 103*0c0fd78eSBarry Smith */ 1048fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 1058fe8eb27SJed Brown { 1068fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1078fe8eb27SJed Brown PetscErrorCode ierr; 1088fe8eb27SJed Brown 1098fe8eb27SJed Brown PetscFunctionBegin; 1108fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);} 1118fe8eb27SJed Brown PetscFunctionReturn(0); 1128fe8eb27SJed Brown } 1138fe8eb27SJed Brown 114*0c0fd78eSBarry Smith /* 115*0c0fd78eSBarry Smith x = diag(right)*x 116*0c0fd78eSBarry Smith */ 1178fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 1188fe8eb27SJed Brown { 1198fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1208fe8eb27SJed Brown PetscErrorCode ierr; 1218fe8eb27SJed Brown 1228fe8eb27SJed Brown PetscFunctionBegin; 1238fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);} 1248fe8eb27SJed Brown PetscFunctionReturn(0); 1258fe8eb27SJed Brown } 1268fe8eb27SJed Brown 127*0c0fd78eSBarry Smith /* 128*0c0fd78eSBarry Smith Y = vscale*Y + diag(dshift)*X + vshift*X 129*0c0fd78eSBarry Smith 130*0c0fd78eSBarry Smith On input Y already contains A*x 131*0c0fd78eSBarry Smith */ 1328fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 1338fe8eb27SJed Brown { 1348fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1358fe8eb27SJed Brown PetscErrorCode ierr; 1368fe8eb27SJed Brown 1378fe8eb27SJed Brown PetscFunctionBegin; 1388fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 1398fe8eb27SJed Brown PetscInt i,m; 1408fe8eb27SJed Brown const PetscScalar *x,*d; 1418fe8eb27SJed Brown PetscScalar *y; 1428fe8eb27SJed Brown ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr); 1438fe8eb27SJed Brown ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1448fe8eb27SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 1458fe8eb27SJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1468fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 1478fe8eb27SJed Brown ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1488fe8eb27SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 1498fe8eb27SJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 150*0c0fd78eSBarry Smith } else { 151027c5fb5SJed Brown ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr); 1528fe8eb27SJed Brown } 153*0c0fd78eSBarry Smith if (shell->vshift) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */ 1548fe8eb27SJed Brown PetscFunctionReturn(0); 1558fe8eb27SJed Brown } 156e51e0e81SBarry Smith 1579d225801SJed Brown /*@ 158a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 159b4fd4287SBarry Smith 160c7fcc2eaSBarry Smith Not Collective 161c7fcc2eaSBarry Smith 162b4fd4287SBarry Smith Input Parameter: 163b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 164b4fd4287SBarry Smith 165b4fd4287SBarry Smith Output Parameter: 166b4fd4287SBarry Smith . ctx - the user provided context 167b4fd4287SBarry Smith 16815091d37SBarry Smith Level: advanced 16915091d37SBarry Smith 170daf670e6SBarry Smith Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this 171daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 172a62d957aSLois Curfman McInnes 173a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context 174a62d957aSLois Curfman McInnes 175ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 1760bc0a719SBarry Smith @*/ 1778fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 178b4fd4287SBarry Smith { 179dfbe8321SBarry Smith PetscErrorCode ierr; 180ace3abfcSBarry Smith PetscBool flg; 181273d9f13SBarry Smith 1823a40ed3dSBarry Smith PetscFunctionBegin; 1830700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1844482741eSBarry Smith PetscValidPointer(ctx,2); 185251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 186940183f3SBarry Smith if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 187ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix"); 1883a40ed3dSBarry Smith PetscFunctionReturn(0); 189b4fd4287SBarry Smith } 190b4fd4287SBarry Smith 191dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 192e51e0e81SBarry Smith { 193dfbe8321SBarry Smith PetscErrorCode ierr; 194bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 195ed3cc1f0SBarry Smith 1963a40ed3dSBarry Smith PetscFunctionBegin; 19728f357e3SAlex Fikl if (shell->ops->destroy) { 19828f357e3SAlex Fikl ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr); 199bf0cc555SLisandro Dalcin } 200*0c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 201*0c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 202*0c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 2038fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 2048fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 2055edf6516SJed Brown ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr); 2065edf6516SJed Brown ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr); 207bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 2083a40ed3dSBarry Smith PetscFunctionReturn(0); 209e51e0e81SBarry Smith } 210e51e0e81SBarry Smith 2117fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str) 2127fabda3fSAlex Fikl { 21328f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data; 2147fabda3fSAlex Fikl PetscErrorCode ierr; 215cb8c8a77SBarry Smith PetscBool matflg; 2167fabda3fSAlex Fikl 2177fabda3fSAlex Fikl PetscFunctionBegin; 21828f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr); 219cb8c8a77SBarry Smith if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL"); 2207fabda3fSAlex Fikl 221cb8c8a77SBarry Smith ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 222cb8c8a77SBarry Smith ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr); 2237fabda3fSAlex Fikl 224cb8c8a77SBarry Smith if (shellA->ops->copy) { 22528f357e3SAlex Fikl ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr); 226cb8c8a77SBarry Smith } 2277fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 2287fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 229*0c0fd78eSBarry Smith if (shellA->dshift) { 230*0c0fd78eSBarry Smith if (!shellB->dshift) { 231*0c0fd78eSBarry Smith ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr); 2327fabda3fSAlex Fikl } 233*0c0fd78eSBarry Smith ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr); 234*0c0fd78eSBarry Smith shellB->dshift = shellB->dshift; 2357fabda3fSAlex Fikl } else { 2367fabda3fSAlex Fikl shellB->dshift = NULL; 2377fabda3fSAlex Fikl } 238*0c0fd78eSBarry Smith if (shellA->left) { 239*0c0fd78eSBarry Smith if (!shellB->left) { 240*0c0fd78eSBarry Smith ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr); 2417fabda3fSAlex Fikl } 242*0c0fd78eSBarry Smith ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr); 243*0c0fd78eSBarry Smith shellB->left = shellB->left; 2447fabda3fSAlex Fikl } else { 2457fabda3fSAlex Fikl shellB->left = NULL; 2467fabda3fSAlex Fikl } 247*0c0fd78eSBarry Smith if (shellA->right) { 248*0c0fd78eSBarry Smith if (!shellB->right) { 249*0c0fd78eSBarry Smith ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr); 2507fabda3fSAlex Fikl } 251*0c0fd78eSBarry Smith ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr); 252*0c0fd78eSBarry Smith shellB->right = shellB->right; 2537fabda3fSAlex Fikl } else { 2547fabda3fSAlex Fikl shellB->right = NULL; 2557fabda3fSAlex Fikl } 2567fabda3fSAlex Fikl PetscFunctionReturn(0); 2577fabda3fSAlex Fikl } 2587fabda3fSAlex Fikl 259cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M) 260cb8c8a77SBarry Smith { 261cb8c8a77SBarry Smith PetscErrorCode ierr; 262cb8c8a77SBarry Smith void *ctx; 263cb8c8a77SBarry Smith 264cb8c8a77SBarry Smith PetscFunctionBegin; 265cb8c8a77SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 266cb8c8a77SBarry Smith ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr); 267cb8c8a77SBarry Smith ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 268cb8c8a77SBarry Smith PetscFunctionReturn(0); 269cb8c8a77SBarry Smith } 270cb8c8a77SBarry Smith 271dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 272ef66eb69SBarry Smith { 273ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 274dfbe8321SBarry Smith PetscErrorCode ierr; 27525578ef6SJed Brown Vec xx; 276e3f487b0SBarry Smith PetscObjectState instate,outstate; 277ef66eb69SBarry Smith 278ef66eb69SBarry Smith PetscFunctionBegin; 2798fe8eb27SJed Brown ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr); 280e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 28128f357e3SAlex Fikl ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr); 282e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 283e3f487b0SBarry Smith if (instate == outstate) { 284e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 285e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 286e3f487b0SBarry Smith } 2878fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2888fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 289ef66eb69SBarry Smith PetscFunctionReturn(0); 290ef66eb69SBarry Smith } 291ef66eb69SBarry Smith 2925edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 2935edf6516SJed Brown { 2945edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2955edf6516SJed Brown PetscErrorCode ierr; 2965edf6516SJed Brown 2975edf6516SJed Brown PetscFunctionBegin; 2985edf6516SJed Brown if (y == z) { 2995edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 3005edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 301b8430d49SBarry Smith ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr); 3025edf6516SJed Brown } else { 3035edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 3045edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 3055edf6516SJed Brown } 3065edf6516SJed Brown PetscFunctionReturn(0); 3075edf6516SJed Brown } 3085edf6516SJed Brown 30918be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 31018be62a5SSatish Balay { 31118be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 31218be62a5SSatish Balay PetscErrorCode ierr; 31325578ef6SJed Brown Vec xx; 314e3f487b0SBarry Smith PetscObjectState instate,outstate; 31518be62a5SSatish Balay 31618be62a5SSatish Balay PetscFunctionBegin; 3178fe8eb27SJed Brown ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr); 318e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 319*0c0fd78eSBarry Smith if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL"); 32028f357e3SAlex Fikl ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr); 321e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 322e3f487b0SBarry Smith if (instate == outstate) { 323e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 324e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 325e3f487b0SBarry Smith } 3268fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 3278fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 32818be62a5SSatish Balay PetscFunctionReturn(0); 32918be62a5SSatish Balay } 33018be62a5SSatish Balay 3315edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 3325edf6516SJed Brown { 3335edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 3345edf6516SJed Brown PetscErrorCode ierr; 3355edf6516SJed Brown 3365edf6516SJed Brown PetscFunctionBegin; 3375edf6516SJed Brown if (y == z) { 3385edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 3395edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 3405edf6516SJed Brown ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr); 3415edf6516SJed Brown } else { 3425edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 3435edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 3445edf6516SJed Brown } 3455edf6516SJed Brown PetscFunctionReturn(0); 3465edf6516SJed Brown } 3475edf6516SJed Brown 348*0c0fd78eSBarry Smith /* 349*0c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 350*0c0fd78eSBarry Smith */ 35118be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 35218be62a5SSatish Balay { 35318be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 35418be62a5SSatish Balay PetscErrorCode ierr; 35518be62a5SSatish Balay 35618be62a5SSatish Balay PetscFunctionBegin; 357*0c0fd78eSBarry Smith if (shell->ops->getdiagonal) { 35828f357e3SAlex Fikl ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr); 359*0c0fd78eSBarry Smith } else { 360*0c0fd78eSBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 361*0c0fd78eSBarry Smith } 36218be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 3638fe8eb27SJed Brown if (shell->dshift) { 364*0c0fd78eSBarry Smith ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr); 36518be62a5SSatish Balay } 366*0c0fd78eSBarry Smith ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 3678fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 3688fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 36918be62a5SSatish Balay PetscFunctionReturn(0); 37018be62a5SSatish Balay } 37118be62a5SSatish Balay 372f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 373ef66eb69SBarry Smith { 374ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3758fe8eb27SJed Brown PetscErrorCode ierr; 376b24ad042SBarry Smith 377ef66eb69SBarry Smith PetscFunctionBegin; 378*0c0fd78eSBarry Smith if (shell->left || shell->right) { 3798fe8eb27SJed Brown if (!shell->dshift) { 380*0c0fd78eSBarry Smith ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr); 381*0c0fd78eSBarry Smith ierr = VecSet(shell->dshift,a);CHKERRQ(ierr); 382*0c0fd78eSBarry Smith } else { 383*0c0fd78eSBarry Smith if (shell->left) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 384*0c0fd78eSBarry Smith if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 385*0c0fd78eSBarry Smith ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 386*0c0fd78eSBarry Smith } 3878fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 3888fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 3898fe8eb27SJed Brown } else shell->vshift += a; 390ef66eb69SBarry Smith PetscFunctionReturn(0); 391ef66eb69SBarry Smith } 392ef66eb69SBarry Smith 3936464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 3946464bf51SAlex Fikl { 3956464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 3966464bf51SAlex Fikl PetscErrorCode ierr; 3976464bf51SAlex Fikl 3986464bf51SAlex Fikl PetscFunctionBegin; 399*0c0fd78eSBarry Smith if (ins == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported with INSERT_VALUES"); 400*0c0fd78eSBarry Smith if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);} 401*0c0fd78eSBarry Smith if (shell->left || shell->right) { 402*0c0fd78eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not yet implemented"); 403*0c0fd78eSBarry Smith } else { 404*0c0fd78eSBarry Smith ierr = VecAXPY(shell->dshift,1.0,D);CHKERRQ(ierr); 4056464bf51SAlex Fikl } 4066464bf51SAlex Fikl PetscFunctionReturn(0); 4076464bf51SAlex Fikl } 4086464bf51SAlex Fikl 409f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 410ef66eb69SBarry Smith { 411ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 4128fe8eb27SJed Brown PetscErrorCode ierr; 413b24ad042SBarry Smith 414ef66eb69SBarry Smith PetscFunctionBegin; 415f4df32b1SMatthew Knepley shell->vscale *= a; 416*0c0fd78eSBarry Smith shell->vshift *= a; 4172205254eSKarl Rupp if (shell->dshift) { 4182205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 419*0c0fd78eSBarry Smith } 4208fe8eb27SJed Brown PetscFunctionReturn(0); 42118be62a5SSatish Balay } 4228fe8eb27SJed Brown 4238fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 4248fe8eb27SJed Brown { 4258fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 4268fe8eb27SJed Brown PetscErrorCode ierr; 4278fe8eb27SJed Brown 4288fe8eb27SJed Brown PetscFunctionBegin; 4298fe8eb27SJed Brown if (left) { 430*0c0fd78eSBarry Smith if (!shell->left) { 431*0c0fd78eSBarry Smith ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr); 4328fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 433*0c0fd78eSBarry Smith } else { 434*0c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 43518be62a5SSatish Balay } 436ef66eb69SBarry Smith } 4378fe8eb27SJed Brown if (right) { 438*0c0fd78eSBarry Smith if (!shell->right) { 439*0c0fd78eSBarry Smith ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr); 4408fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 441*0c0fd78eSBarry Smith } else { 442*0c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 4438fe8eb27SJed Brown } 4448fe8eb27SJed Brown } 445ef66eb69SBarry Smith PetscFunctionReturn(0); 446ef66eb69SBarry Smith } 447ef66eb69SBarry Smith 448dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 449ef66eb69SBarry Smith { 450ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 451*0c0fd78eSBarry Smith PetscErrorCode ierr; 452ef66eb69SBarry Smith 453ef66eb69SBarry Smith PetscFunctionBegin; 4548fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 455ef66eb69SBarry Smith shell->vshift = 0.0; 456ef66eb69SBarry Smith shell->vscale = 1.0; 457*0c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 458*0c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 459*0c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 460ef66eb69SBarry Smith } 461ef66eb69SBarry Smith PetscFunctionReturn(0); 462ef66eb69SBarry Smith } 463ef66eb69SBarry Smith 464cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*); 465b951964fSBarry Smith 4663b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 4673b49f96aSBarry Smith { 4683b49f96aSBarry Smith PetscFunctionBegin; 4693b49f96aSBarry Smith *missing = PETSC_FALSE; 4703b49f96aSBarry Smith PetscFunctionReturn(0); 4713b49f96aSBarry Smith } 4723b49f96aSBarry Smith 47309dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 47420563c6bSBarry Smith 0, 47520563c6bSBarry Smith 0, 476*0c0fd78eSBarry Smith MatMult_Shell, 477*0c0fd78eSBarry Smith /* 4*/ MatMultAdd_Shell, 478*0c0fd78eSBarry Smith MatMultTranspose_Shell, 479*0c0fd78eSBarry Smith MatMultTransposeAdd_Shell, 480b951964fSBarry Smith 0, 481b951964fSBarry Smith 0, 482b951964fSBarry Smith 0, 48397304618SKris Buschelman /*10*/ 0, 484b951964fSBarry Smith 0, 485b951964fSBarry Smith 0, 486b951964fSBarry Smith 0, 487b951964fSBarry Smith 0, 48897304618SKris Buschelman /*15*/ 0, 489b951964fSBarry Smith 0, 490*0c0fd78eSBarry Smith MatGetDiagonal_Shell, 4918fe8eb27SJed Brown MatDiagonalScale_Shell, 492b951964fSBarry Smith 0, 49397304618SKris Buschelman /*20*/ 0, 494ef66eb69SBarry Smith MatAssemblyEnd_Shell, 495b951964fSBarry Smith 0, 496b951964fSBarry Smith 0, 497d519adbfSMatthew Knepley /*24*/ 0, 498b951964fSBarry Smith 0, 499b951964fSBarry Smith 0, 500b951964fSBarry Smith 0, 501b951964fSBarry Smith 0, 502d519adbfSMatthew Knepley /*29*/ 0, 503b951964fSBarry Smith 0, 504273d9f13SBarry Smith 0, 505b951964fSBarry Smith 0, 506b951964fSBarry Smith 0, 507cb8c8a77SBarry Smith /*34*/ MatDuplicate_Shell, 508b951964fSBarry Smith 0, 509b951964fSBarry Smith 0, 51009dc0095SBarry Smith 0, 51109dc0095SBarry Smith 0, 512d519adbfSMatthew Knepley /*39*/ 0, 51309dc0095SBarry Smith 0, 51409dc0095SBarry Smith 0, 51509dc0095SBarry Smith 0, 516cb8c8a77SBarry Smith MatCopy_Shell, 517d519adbfSMatthew Knepley /*44*/ 0, 518ef66eb69SBarry Smith MatScale_Shell, 519ef66eb69SBarry Smith MatShift_Shell, 5206464bf51SAlex Fikl MatDiagonalSet_Shell, 52109dc0095SBarry Smith 0, 522f73d5cc4SBarry Smith /*49*/ 0, 52309dc0095SBarry Smith 0, 52409dc0095SBarry Smith 0, 52509dc0095SBarry Smith 0, 52609dc0095SBarry Smith 0, 527d519adbfSMatthew Knepley /*54*/ 0, 52809dc0095SBarry Smith 0, 52909dc0095SBarry Smith 0, 53009dc0095SBarry Smith 0, 53109dc0095SBarry Smith 0, 532d519adbfSMatthew Knepley /*59*/ 0, 533b9b97703SBarry Smith MatDestroy_Shell, 53409dc0095SBarry Smith 0, 535357abbc8SBarry Smith 0, 536273d9f13SBarry Smith 0, 537d519adbfSMatthew Knepley /*64*/ 0, 538273d9f13SBarry Smith 0, 539273d9f13SBarry Smith 0, 540273d9f13SBarry Smith 0, 541273d9f13SBarry Smith 0, 542d519adbfSMatthew Knepley /*69*/ 0, 543273d9f13SBarry Smith 0, 544c87e5d42SMatthew Knepley MatConvert_Shell, 545273d9f13SBarry Smith 0, 54697304618SKris Buschelman 0, 547d519adbfSMatthew Knepley /*74*/ 0, 54897304618SKris Buschelman 0, 54997304618SKris Buschelman 0, 55097304618SKris Buschelman 0, 55197304618SKris Buschelman 0, 552d519adbfSMatthew Knepley /*79*/ 0, 55397304618SKris Buschelman 0, 55497304618SKris Buschelman 0, 55597304618SKris Buschelman 0, 556865e5f61SKris Buschelman 0, 557d519adbfSMatthew Knepley /*84*/ 0, 558865e5f61SKris Buschelman 0, 559865e5f61SKris Buschelman 0, 560865e5f61SKris Buschelman 0, 561865e5f61SKris Buschelman 0, 562d519adbfSMatthew Knepley /*89*/ 0, 563865e5f61SKris Buschelman 0, 564865e5f61SKris Buschelman 0, 565865e5f61SKris Buschelman 0, 566865e5f61SKris Buschelman 0, 567d519adbfSMatthew Knepley /*94*/ 0, 568865e5f61SKris Buschelman 0, 569865e5f61SKris Buschelman 0, 5703964eb88SJed Brown 0, 5713964eb88SJed Brown 0, 5723964eb88SJed Brown /*99*/ 0, 5733964eb88SJed Brown 0, 5743964eb88SJed Brown 0, 5753964eb88SJed Brown 0, 5763964eb88SJed Brown 0, 5773964eb88SJed Brown /*104*/ 0, 5783964eb88SJed Brown 0, 5793964eb88SJed Brown 0, 5803964eb88SJed Brown 0, 5813964eb88SJed Brown 0, 5823964eb88SJed Brown /*109*/ 0, 5833964eb88SJed Brown 0, 5843964eb88SJed Brown 0, 5853964eb88SJed Brown 0, 5863b49f96aSBarry Smith MatMissingDiagonal_Shell, 5873964eb88SJed Brown /*114*/ 0, 5883964eb88SJed Brown 0, 5893964eb88SJed Brown 0, 5903964eb88SJed Brown 0, 5913964eb88SJed Brown 0, 5923964eb88SJed Brown /*119*/ 0, 5933964eb88SJed Brown 0, 5943964eb88SJed Brown 0, 5953964eb88SJed Brown 0, 5963964eb88SJed Brown 0, 5973964eb88SJed Brown /*124*/ 0, 5983964eb88SJed Brown 0, 5993964eb88SJed Brown 0, 6003964eb88SJed Brown 0, 6013964eb88SJed Brown 0, 6023964eb88SJed Brown /*129*/ 0, 6033964eb88SJed Brown 0, 6043964eb88SJed Brown 0, 6053964eb88SJed Brown 0, 6063964eb88SJed Brown 0, 6073964eb88SJed Brown /*134*/ 0, 6083964eb88SJed Brown 0, 6093964eb88SJed Brown 0, 6103964eb88SJed Brown 0, 6113964eb88SJed Brown 0, 6123964eb88SJed Brown /*139*/ 0, 613f9426fe0SMark Adams 0, 6143964eb88SJed Brown 0 6153964eb88SJed Brown }; 616273d9f13SBarry Smith 6170bad9183SKris Buschelman /*MC 618fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 6190bad9183SKris Buschelman 6200bad9183SKris Buschelman Level: advanced 6210bad9183SKris Buschelman 622*0c0fd78eSBarry Smith .seealso: MatCreateShell() 6230bad9183SKris Buschelman M*/ 6240bad9183SKris Buschelman 6258cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 626273d9f13SBarry Smith { 627273d9f13SBarry Smith Mat_Shell *b; 628dfbe8321SBarry Smith PetscErrorCode ierr; 629273d9f13SBarry Smith 630273d9f13SBarry Smith PetscFunctionBegin; 631273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 632273d9f13SBarry Smith 633b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 634273d9f13SBarry Smith A->data = (void*)b; 635273d9f13SBarry Smith 63626283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 63726283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 638273d9f13SBarry Smith 639273d9f13SBarry Smith b->ctx = 0; 640ef66eb69SBarry Smith b->vshift = 0.0; 641ef66eb69SBarry Smith b->vscale = 1.0; 642*0c0fd78eSBarry Smith b->managescalingshifts = PETSC_TRUE; 643273d9f13SBarry Smith A->assembled = PETSC_TRUE; 644210f0121SBarry Smith A->preallocated = PETSC_FALSE; 6452205254eSKarl Rupp 64617667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 647273d9f13SBarry Smith PetscFunctionReturn(0); 648273d9f13SBarry Smith } 649e51e0e81SBarry Smith 6504b828684SBarry Smith /*@C 651052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 652ff756334SLois Curfman McInnes private data storage format. 653e51e0e81SBarry Smith 654c7fcc2eaSBarry Smith Collective on MPI_Comm 655c7fcc2eaSBarry Smith 656e51e0e81SBarry Smith Input Parameters: 657c7fcc2eaSBarry Smith + comm - MPI communicator 658c7fcc2eaSBarry Smith . m - number of local rows (must be given) 659c7fcc2eaSBarry Smith . n - number of local columns (must be given) 660c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 661c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 662c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 663e51e0e81SBarry Smith 664ff756334SLois Curfman McInnes Output Parameter: 66544cd7ae7SLois Curfman McInnes . A - the matrix 666e51e0e81SBarry Smith 667ff2fd236SBarry Smith Level: advanced 668ff2fd236SBarry Smith 669f39d1f56SLois Curfman McInnes Usage: 6707b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 671f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 672c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 673f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 674f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 675f39d1f56SLois Curfman McInnes 676ff756334SLois Curfman McInnes Notes: 677ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 678ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 679ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 680e51e0e81SBarry Smith 681daf670e6SBarry Smith Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this 682daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 683daf670e6SBarry Smith in as the ctx argument. 684f38a8d56SBarry Smith 685f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 686f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 687645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 688645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 689645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 690645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 691645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 692645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 693645985a0SLois Curfman McInnes For example, 694f39d1f56SLois Curfman McInnes 695*0c0fd78eSBarry Smith MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), and MatScale() internally so these 696*0c0fd78eSBarry Smith operations cannot be overwritten unless MatShellSetManageScalingShifts() is called. 697*0c0fd78eSBarry Smith 698f39d1f56SLois Curfman McInnes $ 699f39d1f56SLois Curfman McInnes $ Vec x, y 7007b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 701645985a0SLois Curfman McInnes $ Mat A 702f39d1f56SLois Curfman McInnes $ 703c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 704c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 705f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 706c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 707c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 708c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 709645985a0SLois Curfman McInnes $ MatMult(A,x,y); 710645985a0SLois Curfman McInnes $ MatDestroy(A); 711f39d1f56SLois Curfman McInnes $ VecDestroy(y); VecDestroy(x); 712645985a0SLois Curfman McInnes $ 713e51e0e81SBarry Smith 714*0c0fd78eSBarry Smith For rectangular matrices do all the scalings and shifts make sense? 715*0c0fd78eSBarry Smith 716*0c0fd78eSBarry Smith Developers Notes: Regarding shifting and scaling. The general form is 717*0c0fd78eSBarry Smith 718*0c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 719*0c0fd78eSBarry Smith 720*0c0fd78eSBarry Smith The order you apply the operations is important. For example if you have a dshift then 721*0c0fd78eSBarry Smith apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift 722*0c0fd78eSBarry Smith you get s*vscale*A + diag(shift) 723*0c0fd78eSBarry Smith 724*0c0fd78eSBarry Smith A is the user provided function. 725*0c0fd78eSBarry Smith 7260b627109SLois Curfman McInnes .keywords: matrix, shell, create 7270b627109SLois Curfman McInnes 728*0c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts() 729e51e0e81SBarry Smith @*/ 7307087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 731e51e0e81SBarry Smith { 732dfbe8321SBarry Smith PetscErrorCode ierr; 733ed3cc1f0SBarry Smith 7343a40ed3dSBarry Smith PetscFunctionBegin; 735f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 736f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 737273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 738273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 7390fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 740273d9f13SBarry Smith PetscFunctionReturn(0); 741c7fcc2eaSBarry Smith } 742c7fcc2eaSBarry Smith 743c6866cfdSSatish Balay /*@ 744273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 745c7fcc2eaSBarry Smith 7463f9fe445SBarry Smith Logically Collective on Mat 747c7fcc2eaSBarry Smith 748273d9f13SBarry Smith Input Parameters: 749273d9f13SBarry Smith + mat - the shell matrix 750273d9f13SBarry Smith - ctx - the context 751273d9f13SBarry Smith 752273d9f13SBarry Smith Level: advanced 753273d9f13SBarry Smith 754daf670e6SBarry Smith Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this 755daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 756273d9f13SBarry Smith 757273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 7580bc0a719SBarry Smith @*/ 7597087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 760273d9f13SBarry Smith { 761273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 762dfbe8321SBarry Smith PetscErrorCode ierr; 763ace3abfcSBarry Smith PetscBool flg; 764273d9f13SBarry Smith 765273d9f13SBarry Smith PetscFunctionBegin; 7660700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 767251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 768273d9f13SBarry Smith if (flg) { 769273d9f13SBarry Smith shell->ctx = ctx; 770ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix"); 7713a40ed3dSBarry Smith PetscFunctionReturn(0); 772e51e0e81SBarry Smith } 773e51e0e81SBarry Smith 774*0c0fd78eSBarry Smith /*@ 775*0c0fd78eSBarry Smith MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately 776*0c0fd78eSBarry Smith after MatCreateShell() 777*0c0fd78eSBarry Smith 778*0c0fd78eSBarry Smith Logically Collective on Mat 779*0c0fd78eSBarry Smith 780*0c0fd78eSBarry Smith Input Parameter: 781*0c0fd78eSBarry Smith . mat - the shell matrix 782*0c0fd78eSBarry Smith 783*0c0fd78eSBarry Smith Level: advanced 784*0c0fd78eSBarry Smith 785*0c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation() 786*0c0fd78eSBarry Smith @*/ 787*0c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A) 788*0c0fd78eSBarry Smith { 789*0c0fd78eSBarry Smith PetscErrorCode ierr; 790*0c0fd78eSBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 791*0c0fd78eSBarry Smith PetscBool flg; 792*0c0fd78eSBarry Smith 793*0c0fd78eSBarry Smith PetscFunctionBegin; 794*0c0fd78eSBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 795*0c0fd78eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr); 796*0c0fd78eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 797*0c0fd78eSBarry Smith shell->managescalingshifts = PETSC_FALSE; 798*0c0fd78eSBarry Smith A->ops->diagonalscale = 0; 799*0c0fd78eSBarry Smith A->ops->scale = 0; 800*0c0fd78eSBarry Smith A->ops->shift = 0; 801*0c0fd78eSBarry Smith A->ops->diagonalset = 0; 802*0c0fd78eSBarry Smith PetscFunctionReturn(0); 803*0c0fd78eSBarry Smith } 804*0c0fd78eSBarry Smith 805c16cb8f2SBarry Smith /*@C 806f3b1f45cSBarry Smith MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function. 807f3b1f45cSBarry Smith 808f3b1f45cSBarry Smith Logically Collective on Mat 809f3b1f45cSBarry Smith 810f3b1f45cSBarry Smith Input Parameters: 811f3b1f45cSBarry Smith + mat - the shell matrix 812f3b1f45cSBarry Smith . f - the function 813f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 814f3b1f45cSBarry Smith - ctx - an optional context for the function 815f3b1f45cSBarry Smith 816f3b1f45cSBarry Smith Output Parameter: 817f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 818f3b1f45cSBarry Smith 819f3b1f45cSBarry Smith Options Database: 820f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 821f3b1f45cSBarry Smith 822f3b1f45cSBarry Smith Level: advanced 823f3b1f45cSBarry Smith 824f3b1f45cSBarry Smith Fortran Notes: Not supported from Fortran 825f3b1f45cSBarry Smith 826f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose() 827f3b1f45cSBarry Smith @*/ 828f3b1f45cSBarry Smith PetscErrorCode MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 829f3b1f45cSBarry Smith { 830f3b1f45cSBarry Smith PetscErrorCode ierr; 831f3b1f45cSBarry Smith PetscInt m,n; 832f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 833f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 83474e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 835f3b1f45cSBarry Smith 836f3b1f45cSBarry Smith PetscFunctionBegin; 837f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 838f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr); 839f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 840f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr); 841f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 842f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 843f3b1f45cSBarry Smith 844f3b1f45cSBarry Smith ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr); 845f3b1f45cSBarry Smith ierr = MatComputeExplicitOperator(mat,&Dmat);CHKERRQ(ierr); 846f3b1f45cSBarry Smith 847f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 848f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 849f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 850f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 851f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 852f3b1f45cSBarry Smith flag = PETSC_FALSE; 853f3b1f45cSBarry Smith if (v) { 854f3b1f45cSBarry 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)); 855f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 856f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 857f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 858f3b1f45cSBarry Smith } 859f3b1f45cSBarry Smith } else if (v) { 860f3b1f45cSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n"); 861f3b1f45cSBarry Smith } 862f3b1f45cSBarry Smith if (flg) *flg = flag; 863f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 864f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 865f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 866f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 867f3b1f45cSBarry Smith PetscFunctionReturn(0); 868f3b1f45cSBarry Smith } 869f3b1f45cSBarry Smith 870f3b1f45cSBarry Smith /*@C 871f3b1f45cSBarry Smith MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function. 872f3b1f45cSBarry Smith 873f3b1f45cSBarry Smith Logically Collective on Mat 874f3b1f45cSBarry Smith 875f3b1f45cSBarry Smith Input Parameters: 876f3b1f45cSBarry Smith + mat - the shell matrix 877f3b1f45cSBarry Smith . f - the function 878f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 879f3b1f45cSBarry Smith - ctx - an optional context for the function 880f3b1f45cSBarry Smith 881f3b1f45cSBarry Smith Output Parameter: 882f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 883f3b1f45cSBarry Smith 884f3b1f45cSBarry Smith Options Database: 885f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 886f3b1f45cSBarry Smith 887f3b1f45cSBarry Smith Level: advanced 888f3b1f45cSBarry Smith 889f3b1f45cSBarry Smith Fortran Notes: Not supported from Fortran 890f3b1f45cSBarry Smith 891f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult() 892f3b1f45cSBarry Smith @*/ 893f3b1f45cSBarry Smith PetscErrorCode MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 894f3b1f45cSBarry Smith { 895f3b1f45cSBarry Smith PetscErrorCode ierr; 896f3b1f45cSBarry Smith Vec x,y,z; 897f3b1f45cSBarry Smith PetscInt m,n,M,N; 898f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 899f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 90074e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 901f3b1f45cSBarry Smith 902f3b1f45cSBarry Smith PetscFunctionBegin; 903f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 904f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr); 905f3b1f45cSBarry Smith ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr); 906f3b1f45cSBarry Smith ierr = VecDuplicate(y,&z);CHKERRQ(ierr); 907f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 908f3b1f45cSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 909f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr); 910f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 911f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 912f3b1f45cSBarry Smith ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr); 913f3b1f45cSBarry Smith ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr); 914f3b1f45cSBarry Smith ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr); 915f3b1f45cSBarry Smith 916f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 917f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 918f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 919f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 920f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 921f3b1f45cSBarry Smith flag = PETSC_FALSE; 922f3b1f45cSBarry Smith if (v) { 923f3b1f45cSBarry 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)); 924f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 925f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 926f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 927f3b1f45cSBarry Smith } 928f3b1f45cSBarry Smith } else if (v) { 929f3b1f45cSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n"); 930f3b1f45cSBarry Smith } 931f3b1f45cSBarry Smith if (flg) *flg = flag; 932f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 933f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 934f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 935f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 936f3b1f45cSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 937f3b1f45cSBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 938f3b1f45cSBarry Smith ierr = VecDestroy(&z);CHKERRQ(ierr); 939f3b1f45cSBarry Smith PetscFunctionReturn(0); 940f3b1f45cSBarry Smith } 941f3b1f45cSBarry Smith 942f3b1f45cSBarry Smith /*@C 9433a3eedf2SBarry Smith MatShellSetOperation - Allows user to set a matrix operation for 9443a3eedf2SBarry Smith a shell matrix. 945*0c0fd78eSBarry Smith MatShellSetOperation - Allows user to set a matrix operation for a shell matrix. 946e51e0e81SBarry Smith 9473f9fe445SBarry Smith Logically Collective on Mat 948fee21e36SBarry Smith 949c7fcc2eaSBarry Smith Input Parameters: 950c7fcc2eaSBarry Smith + mat - the shell matrix 951c7fcc2eaSBarry Smith . op - the name of the operation 952c7fcc2eaSBarry Smith - f - the function that provides the operation. 953c7fcc2eaSBarry Smith 95415091d37SBarry Smith Level: advanced 95515091d37SBarry Smith 956fae171e0SBarry Smith Usage: 957e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 958f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 959c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 9600b627109SLois Curfman McInnes 961a62d957aSLois Curfman McInnes Notes: 962e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 9631c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 964a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 9651c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 966a62d957aSLois Curfman McInnes 967a4d85fd7SBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 968deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 969deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 970deebb3c3SLois Curfman McInnes routines, e.g., 971a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 972a62d957aSLois Curfman McInnes 9734aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 9744aa34b0aSBarry Smith nonzero on failure. 9754aa34b0aSBarry Smith 976a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 977a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 978a62d957aSLois Curfman McInnes set by MatCreateShell(). 979a62d957aSLois Curfman McInnes 9802a7a6963SBarry Smith Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 981501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 982501d9185SBarry Smith 983*0c0fd78eSBarry Smith Use MatSetOperation() to set an operation for any matrix type 984*0c0fd78eSBarry Smith 985a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation 986a62d957aSLois Curfman McInnes 987*0c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts() 988e51e0e81SBarry Smith @*/ 9897087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 990e51e0e81SBarry Smith { 991dfbe8321SBarry Smith PetscErrorCode ierr; 992ace3abfcSBarry Smith PetscBool flg; 993*0c0fd78eSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 994273d9f13SBarry Smith 9953a40ed3dSBarry Smith PetscFunctionBegin; 9960700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 997*0c0fd78eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 998*0c0fd78eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 9995edf6516SJed Brown switch (op) { 10005edf6516SJed Brown case MATOP_DESTROY: 100128f357e3SAlex Fikl shell->ops->destroy = (PetscErrorCode (*)(Mat))f; 10025edf6516SJed Brown break; 10036464bf51SAlex Fikl case MATOP_DIAGONAL_SET: 1004*0c0fd78eSBarry Smith case MATOP_DIAGONAL_SCALE: 1005*0c0fd78eSBarry Smith case MATOP_SHIFT: 1006*0c0fd78eSBarry Smith case MATOP_SCALE: 1007*0c0fd78eSBarry Smith if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()"); 1008*0c0fd78eSBarry Smith (((void(**)(void))mat->ops)[op]) = f; 10096464bf51SAlex Fikl break; 1010*0c0fd78eSBarry Smith case MATOP_GET_DIAGONAL: 1011*0c0fd78eSBarry Smith if (shell->managescalingshifts) shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1012*0c0fd78eSBarry Smith else mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 10135edf6516SJed Brown case MATOP_VIEW: 101440a5a6c4SBarry Smith if (!mat->ops->viewnative) { 101540a5a6c4SBarry Smith mat->ops->viewnative = mat->ops->view; 101640a5a6c4SBarry Smith } 10175edf6516SJed Brown mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 10185edf6516SJed Brown break; 10195edf6516SJed Brown case MATOP_MULT: 1020*0c0fd78eSBarry Smith if (shell->managescalingshifts) shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1021*0c0fd78eSBarry Smith else mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 10225edf6516SJed Brown break; 10235edf6516SJed Brown case MATOP_MULT_TRANSPOSE: 1024*0c0fd78eSBarry Smith if (shell->managescalingshifts) shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1025*0c0fd78eSBarry Smith else mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 10265edf6516SJed Brown break; 10275edf6516SJed Brown default: 10285edf6516SJed Brown (((void(**)(void))mat->ops)[op]) = f; 1029a62d957aSLois Curfman McInnes } 10303a40ed3dSBarry Smith PetscFunctionReturn(0); 1031e51e0e81SBarry Smith } 1032f0479e8cSBarry Smith 1033d4bb536fSBarry Smith /*@C 1034d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 1035d4bb536fSBarry Smith 1036c7fcc2eaSBarry Smith Not Collective 1037c7fcc2eaSBarry Smith 1038d4bb536fSBarry Smith Input Parameters: 1039c7fcc2eaSBarry Smith + mat - the shell matrix 1040c7fcc2eaSBarry Smith - op - the name of the operation 1041d4bb536fSBarry Smith 1042d4bb536fSBarry Smith Output Parameter: 1043d4bb536fSBarry Smith . f - the function that provides the operation. 1044d4bb536fSBarry Smith 104515091d37SBarry Smith Level: advanced 104615091d37SBarry Smith 1047d4bb536fSBarry Smith Notes: 1048e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 1049d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 1050d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 1051d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 1052d4bb536fSBarry Smith 1053d4bb536fSBarry Smith All user-provided functions have the same calling 1054d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 1055d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 1056d4bb536fSBarry Smith routines, e.g., 1057d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1058d4bb536fSBarry Smith 1059d4bb536fSBarry Smith Within each user-defined routine, the user should call 1060d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 1061d4bb536fSBarry Smith set by MatCreateShell(). 1062d4bb536fSBarry Smith 1063d4bb536fSBarry Smith .keywords: matrix, shell, set, operation 1064d4bb536fSBarry Smith 1065ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 1066d4bb536fSBarry Smith @*/ 10677087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 1068d4bb536fSBarry Smith { 1069dfbe8321SBarry Smith PetscErrorCode ierr; 1070ace3abfcSBarry Smith PetscBool flg; 1071273d9f13SBarry Smith 10723a40ed3dSBarry Smith PetscFunctionBegin; 10730700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 107428f357e3SAlex Fikl switch (op) { 107528f357e3SAlex Fikl case MATOP_DESTROY: 1076251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 1077273d9f13SBarry Smith if (flg) { 1078d4bb536fSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 107928f357e3SAlex Fikl *f = (void (*)(void))shell->ops->destroy; 1080c7fcc2eaSBarry Smith } else { 1081c134de8dSSatish Balay *f = (void (*)(void))mat->ops->destroy; 1082d4bb536fSBarry Smith } 108328f357e3SAlex Fikl break; 108428f357e3SAlex Fikl case MATOP_DIAGONAL_SET: 108528f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 108628f357e3SAlex Fikl if (flg) { 108728f357e3SAlex Fikl Mat_Shell *shell = (Mat_Shell*)mat->data; 108828f357e3SAlex Fikl *f = (void (*)(void))shell->ops->diagonalset; 1089c7fcc2eaSBarry Smith } else { 109028f357e3SAlex Fikl *f = (void (*)(void))mat->ops->diagonalset; 109128f357e3SAlex Fikl } 109228f357e3SAlex Fikl break; 109328f357e3SAlex Fikl case MATOP_VIEW: 109428f357e3SAlex Fikl *f = (void (*)(void))mat->ops->view; 109528f357e3SAlex Fikl break; 109628f357e3SAlex Fikl default: 1097c134de8dSSatish Balay *f = (((void (**)(void))mat->ops)[op]); 1098d4bb536fSBarry Smith } 10993a40ed3dSBarry Smith PetscFunctionReturn(0); 1100d4bb536fSBarry Smith } 1101d4bb536fSBarry Smith 1102*0c0fd78eSBarry Smith /*@C 1103*0c0fd78eSBarry Smith MatSetOperation - Allows user to set a matrix operation for any matrix type 1104*0c0fd78eSBarry Smith 1105*0c0fd78eSBarry Smith Logically Collective on Mat 1106*0c0fd78eSBarry Smith 1107*0c0fd78eSBarry Smith Input Parameters: 1108*0c0fd78eSBarry Smith + mat - the shell matrix 1109*0c0fd78eSBarry Smith . op - the name of the operation 1110*0c0fd78eSBarry Smith - f - the function that provides the operation. 1111*0c0fd78eSBarry Smith 1112*0c0fd78eSBarry Smith Level: developer 1113*0c0fd78eSBarry Smith 1114*0c0fd78eSBarry Smith Usage: 1115*0c0fd78eSBarry Smith $ extern PetscErrorCode usermult(Mat,Vec,Vec); 1116*0c0fd78eSBarry Smith $ ierr = MatCreateXXX(comm,...&A); 1117*0c0fd78eSBarry Smith $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 1118*0c0fd78eSBarry Smith 1119*0c0fd78eSBarry Smith Notes: 1120*0c0fd78eSBarry Smith See the file include/petscmat.h for a complete list of matrix 1121*0c0fd78eSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 1122*0c0fd78eSBarry Smith <OPERATION> is the name (in all capital letters) of the 1123*0c0fd78eSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 1124*0c0fd78eSBarry Smith 1125*0c0fd78eSBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 1126*0c0fd78eSBarry Smith sequence as the usual matrix interface routines, since they 1127*0c0fd78eSBarry Smith are intended to be accessed via the usual matrix interface 1128*0c0fd78eSBarry Smith routines, e.g., 1129*0c0fd78eSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1130*0c0fd78eSBarry Smith 1131*0c0fd78eSBarry Smith In particular each function MUST return an error code of 0 on success and 1132*0c0fd78eSBarry Smith nonzero on failure. 1133*0c0fd78eSBarry Smith 1134*0c0fd78eSBarry Smith This routine is distinct from MatShellSetOperation() in that it can be called on any matrix type. 1135*0c0fd78eSBarry Smith 1136*0c0fd78eSBarry Smith .keywords: matrix, shell, set, operation 1137*0c0fd78eSBarry Smith 1138*0c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation() 1139*0c0fd78eSBarry Smith @*/ 1140*0c0fd78eSBarry Smith PetscErrorCode MatSetOperation(Mat mat,MatOperation op,void (*f)(void)) 1141*0c0fd78eSBarry Smith { 1142*0c0fd78eSBarry Smith PetscFunctionBegin; 1143*0c0fd78eSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1144*0c0fd78eSBarry Smith (((void(**)(void))mat->ops)[op]) = f; 1145*0c0fd78eSBarry Smith PetscFunctionReturn(0); 1146*0c0fd78eSBarry Smith } 1147