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 8b45d2f2cSJed Brown #include <petsc-private/matimpl.h> /*I "petscmat.h" I*/ 9b45d2f2cSJed Brown #include <petsc-private/vecimpl.h> 10e51e0e81SBarry Smith 1120563c6bSBarry Smith typedef struct { 126849ba73SBarry Smith PetscErrorCode (*destroy)(Mat); 136849ba73SBarry Smith PetscErrorCode (*mult)(Mat,Vec,Vec); 1418be62a5SSatish Balay PetscErrorCode (*multtranspose)(Mat,Vec,Vec); 1518be62a5SSatish Balay PetscErrorCode (*getdiagonal)(Mat,Vec); 162205254eSKarl Rupp 17ef66eb69SBarry Smith PetscScalar vscale,vshift; 188fe8eb27SJed Brown Vec dshift; 198fe8eb27SJed Brown Vec left,right; 208fe8eb27SJed Brown Vec dshift_owned,left_owned,right_owned; 218fe8eb27SJed Brown Vec left_work,right_work; 225edf6516SJed Brown Vec left_add_work,right_add_work; 238fe8eb27SJed Brown PetscBool usingscaled; 2420563c6bSBarry Smith void *ctx; 2588cf3e7dSBarry Smith } Mat_Shell; 268fe8eb27SJed Brown /* 278fe8eb27SJed Brown The most general expression for the matrix is 288fe8eb27SJed Brown 298fe8eb27SJed Brown S = L (a A + B) R 308fe8eb27SJed Brown 318fe8eb27SJed Brown where 328fe8eb27SJed Brown A is the matrix defined by the user's function 338fe8eb27SJed Brown a is a scalar multiple 348fe8eb27SJed Brown L is left scaling 358fe8eb27SJed Brown R is right scaling 368fe8eb27SJed Brown B is a diagonal shift defined by 378fe8eb27SJed Brown diag(dshift) if the vector dshift is non-NULL 388fe8eb27SJed Brown vscale*identity otherwise 398fe8eb27SJed Brown 408fe8eb27SJed Brown The following identities apply: 418fe8eb27SJed Brown 428fe8eb27SJed Brown Scale by c: 438fe8eb27SJed Brown c [L (a A + B) R] = L [(a c) A + c B] R 448fe8eb27SJed Brown 458fe8eb27SJed Brown Shift by c: 468fe8eb27SJed Brown [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R 478fe8eb27SJed Brown 488fe8eb27SJed Brown Diagonal scaling is achieved by simply multiplying with existing L and R vectors 498fe8eb27SJed Brown 508fe8eb27SJed Brown In the data structure: 518fe8eb27SJed Brown 528fe8eb27SJed Brown vscale=1.0 means no special scaling will be applied 538fe8eb27SJed Brown dshift=NULL means a constant diagonal shift (fall back to vshift) 548fe8eb27SJed Brown vshift=0.0 means no constant diagonal shift, note that vshift is only used if dshift is NULL 558fe8eb27SJed Brown */ 568fe8eb27SJed Brown 578fe8eb27SJed Brown static PetscErrorCode MatMult_Shell(Mat,Vec,Vec); 588fe8eb27SJed Brown static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec); 598fe8eb27SJed Brown static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec); 608fe8eb27SJed Brown 618fe8eb27SJed Brown #undef __FUNCT__ 628fe8eb27SJed Brown #define __FUNCT__ "MatShellUseScaledMethods" 638fe8eb27SJed Brown static PetscErrorCode MatShellUseScaledMethods(Mat Y) 648fe8eb27SJed Brown { 658fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 668fe8eb27SJed Brown 678fe8eb27SJed Brown PetscFunctionBegin; 688fe8eb27SJed Brown if (shell->usingscaled) PetscFunctionReturn(0); 698fe8eb27SJed Brown shell->mult = Y->ops->mult; 708fe8eb27SJed Brown Y->ops->mult = MatMult_Shell; 718fe8eb27SJed Brown if (Y->ops->multtranspose) { 728fe8eb27SJed Brown shell->multtranspose = Y->ops->multtranspose; 738fe8eb27SJed Brown Y->ops->multtranspose = MatMultTranspose_Shell; 748fe8eb27SJed Brown } 758fe8eb27SJed Brown if (Y->ops->getdiagonal) { 768fe8eb27SJed Brown shell->getdiagonal = Y->ops->getdiagonal; 778fe8eb27SJed Brown Y->ops->getdiagonal = MatGetDiagonal_Shell; 788fe8eb27SJed Brown } 798fe8eb27SJed Brown shell->usingscaled = PETSC_TRUE; 808fe8eb27SJed Brown PetscFunctionReturn(0); 818fe8eb27SJed Brown } 828fe8eb27SJed Brown 838fe8eb27SJed Brown #undef __FUNCT__ 848fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleLeft" 858fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 868fe8eb27SJed Brown { 878fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 888fe8eb27SJed Brown PetscErrorCode ierr; 898fe8eb27SJed Brown 908fe8eb27SJed Brown PetscFunctionBegin; 910298fd71SBarry Smith *xx = NULL; 928fe8eb27SJed Brown if (!shell->left) { 938fe8eb27SJed Brown *xx = x; 948fe8eb27SJed Brown } else { 958fe8eb27SJed Brown if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);} 968fe8eb27SJed Brown ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr); 978fe8eb27SJed Brown *xx = shell->left_work; 988fe8eb27SJed Brown } 998fe8eb27SJed Brown PetscFunctionReturn(0); 1008fe8eb27SJed Brown } 1018fe8eb27SJed Brown 1028fe8eb27SJed Brown #undef __FUNCT__ 1038fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleRight" 1048fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 1058fe8eb27SJed Brown { 1068fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1078fe8eb27SJed Brown PetscErrorCode ierr; 1088fe8eb27SJed Brown 1098fe8eb27SJed Brown PetscFunctionBegin; 1100298fd71SBarry Smith *xx = NULL; 1118fe8eb27SJed Brown if (!shell->right) { 1128fe8eb27SJed Brown *xx = x; 1138fe8eb27SJed Brown } else { 1148fe8eb27SJed Brown if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);} 1158fe8eb27SJed Brown ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr); 1168fe8eb27SJed Brown *xx = shell->right_work; 1178fe8eb27SJed Brown } 1188fe8eb27SJed Brown PetscFunctionReturn(0); 1198fe8eb27SJed Brown } 1208fe8eb27SJed Brown 1218fe8eb27SJed Brown #undef __FUNCT__ 1228fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleLeft" 1238fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 1248fe8eb27SJed Brown { 1258fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1268fe8eb27SJed Brown PetscErrorCode ierr; 1278fe8eb27SJed Brown 1288fe8eb27SJed Brown PetscFunctionBegin; 1298fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);} 1308fe8eb27SJed Brown PetscFunctionReturn(0); 1318fe8eb27SJed Brown } 1328fe8eb27SJed Brown 1338fe8eb27SJed Brown #undef __FUNCT__ 1348fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleRight" 1358fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 1368fe8eb27SJed Brown { 1378fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1388fe8eb27SJed Brown PetscErrorCode ierr; 1398fe8eb27SJed Brown 1408fe8eb27SJed Brown PetscFunctionBegin; 1418fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);} 1428fe8eb27SJed Brown PetscFunctionReturn(0); 1438fe8eb27SJed Brown } 1448fe8eb27SJed Brown 1458fe8eb27SJed Brown #undef __FUNCT__ 1468fe8eb27SJed Brown #define __FUNCT__ "MatShellShiftAndScale" 1478fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 1488fe8eb27SJed Brown { 1498fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1508fe8eb27SJed Brown PetscErrorCode ierr; 1518fe8eb27SJed Brown 1528fe8eb27SJed Brown PetscFunctionBegin; 1538fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 1548fe8eb27SJed Brown PetscInt i,m; 1558fe8eb27SJed Brown const PetscScalar *x,*d; 1568fe8eb27SJed Brown PetscScalar *y; 1578fe8eb27SJed Brown ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr); 1588fe8eb27SJed Brown ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1598fe8eb27SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 1608fe8eb27SJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1618fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 1628fe8eb27SJed Brown ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1638fe8eb27SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 1648fe8eb27SJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 165880538faSHong Zhang } else if (PetscAbsScalar(shell->vshift) != 0) { 1668fe8eb27SJed Brown ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr); 167027c5fb5SJed Brown } else if (shell->vscale != 1.0) { 168027c5fb5SJed Brown ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr); 1698fe8eb27SJed Brown } 1708fe8eb27SJed Brown PetscFunctionReturn(0); 1718fe8eb27SJed Brown } 172e51e0e81SBarry Smith 173711e205bSSatish Balay #undef __FUNCT__ 174711e205bSSatish Balay #define __FUNCT__ "MatShellGetContext" 1759d225801SJed Brown /*@ 176a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 177b4fd4287SBarry Smith 178c7fcc2eaSBarry Smith Not Collective 179c7fcc2eaSBarry Smith 180b4fd4287SBarry Smith Input Parameter: 181b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 182b4fd4287SBarry Smith 183b4fd4287SBarry Smith Output Parameter: 184b4fd4287SBarry Smith . ctx - the user provided context 185b4fd4287SBarry Smith 18615091d37SBarry Smith Level: advanced 18715091d37SBarry Smith 188a62d957aSLois Curfman McInnes Notes: 189a62d957aSLois Curfman McInnes This routine is intended for use within various shell matrix routines, 190a62d957aSLois Curfman McInnes as set with MatShellSetOperation(). 191a62d957aSLois Curfman McInnes 192a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context 193a62d957aSLois Curfman McInnes 194ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 1950bc0a719SBarry Smith @*/ 1968fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 197b4fd4287SBarry Smith { 198dfbe8321SBarry Smith PetscErrorCode ierr; 199ace3abfcSBarry Smith PetscBool flg; 200273d9f13SBarry Smith 2013a40ed3dSBarry Smith PetscFunctionBegin; 2020700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2034482741eSBarry Smith PetscValidPointer(ctx,2); 204251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 205940183f3SBarry Smith if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 206ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix"); 2073a40ed3dSBarry Smith PetscFunctionReturn(0); 208b4fd4287SBarry Smith } 209b4fd4287SBarry Smith 210711e205bSSatish Balay #undef __FUNCT__ 211711e205bSSatish Balay #define __FUNCT__ "MatDestroy_Shell" 212dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 213e51e0e81SBarry Smith { 214dfbe8321SBarry Smith PetscErrorCode ierr; 215bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 216ed3cc1f0SBarry Smith 2173a40ed3dSBarry Smith PetscFunctionBegin; 218bf0cc555SLisandro Dalcin if (shell->destroy) { 219bf0cc555SLisandro Dalcin ierr = (*shell->destroy)(mat);CHKERRQ(ierr); 220bf0cc555SLisandro Dalcin } 2218fe8eb27SJed Brown ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr); 2228fe8eb27SJed Brown ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr); 2238fe8eb27SJed Brown ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr); 2248fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 2258fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 2265edf6516SJed Brown ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr); 2275edf6516SJed Brown ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr); 228bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 2293a40ed3dSBarry Smith PetscFunctionReturn(0); 230e51e0e81SBarry Smith } 231e51e0e81SBarry Smith 232711e205bSSatish Balay #undef __FUNCT__ 233ef66eb69SBarry Smith #define __FUNCT__ "MatMult_Shell" 234dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 235ef66eb69SBarry Smith { 236ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 237dfbe8321SBarry Smith PetscErrorCode ierr; 23825578ef6SJed Brown Vec xx; 239ef66eb69SBarry Smith 240ef66eb69SBarry Smith PetscFunctionBegin; 2418fe8eb27SJed Brown ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr); 2428fe8eb27SJed Brown ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr); 2438fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2448fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 245ef66eb69SBarry Smith PetscFunctionReturn(0); 246ef66eb69SBarry Smith } 247ef66eb69SBarry Smith 248ef66eb69SBarry Smith #undef __FUNCT__ 2495edf6516SJed Brown #define __FUNCT__ "MatMultAdd_Shell" 2505edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 2515edf6516SJed Brown { 2525edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2535edf6516SJed Brown PetscErrorCode ierr; 2545edf6516SJed Brown 2555edf6516SJed Brown PetscFunctionBegin; 2565edf6516SJed Brown if (y == z) { 2575edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 2585edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 2595edf6516SJed Brown ierr = VecWAXPY(z,1.0,shell->right_add_work,y);CHKERRQ(ierr); 2605edf6516SJed Brown } else { 2615edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 2625edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 2635edf6516SJed Brown } 2645edf6516SJed Brown PetscFunctionReturn(0); 2655edf6516SJed Brown } 2665edf6516SJed Brown 2675edf6516SJed Brown #undef __FUNCT__ 26818be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell" 26918be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 27018be62a5SSatish Balay { 27118be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 27218be62a5SSatish Balay PetscErrorCode ierr; 27325578ef6SJed Brown Vec xx; 27418be62a5SSatish Balay 27518be62a5SSatish Balay PetscFunctionBegin; 2768fe8eb27SJed Brown ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr); 2778fe8eb27SJed Brown ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr); 2788fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2798fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 28018be62a5SSatish Balay PetscFunctionReturn(0); 28118be62a5SSatish Balay } 28218be62a5SSatish Balay 28318be62a5SSatish Balay #undef __FUNCT__ 2845edf6516SJed Brown #define __FUNCT__ "MatMultTransposeAdd_Shell" 2855edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 2865edf6516SJed Brown { 2875edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2885edf6516SJed Brown PetscErrorCode ierr; 2895edf6516SJed Brown 2905edf6516SJed Brown PetscFunctionBegin; 2915edf6516SJed Brown if (y == z) { 2925edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 2935edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 2945edf6516SJed Brown ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr); 2955edf6516SJed Brown } else { 2965edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 2975edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 2985edf6516SJed Brown } 2995edf6516SJed Brown PetscFunctionReturn(0); 3005edf6516SJed Brown } 3015edf6516SJed Brown 3025edf6516SJed Brown #undef __FUNCT__ 30318be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell" 30418be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 30518be62a5SSatish Balay { 30618be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 30718be62a5SSatish Balay PetscErrorCode ierr; 30818be62a5SSatish Balay 30918be62a5SSatish Balay PetscFunctionBegin; 31018be62a5SSatish Balay ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr); 31118be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 3128fe8eb27SJed Brown if (shell->dshift) { 3138fe8eb27SJed Brown ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr); 3148fe8eb27SJed Brown } else { 31518be62a5SSatish Balay ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 31618be62a5SSatish Balay } 3178fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 3188fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 31918be62a5SSatish Balay PetscFunctionReturn(0); 32018be62a5SSatish Balay } 32118be62a5SSatish Balay 32218be62a5SSatish Balay #undef __FUNCT__ 323ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell" 324f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 325ef66eb69SBarry Smith { 326ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3278fe8eb27SJed Brown PetscErrorCode ierr; 328b24ad042SBarry Smith 329ef66eb69SBarry Smith PetscFunctionBegin; 3308fe8eb27SJed Brown if (shell->left || shell->right || shell->dshift) { 3318fe8eb27SJed Brown if (!shell->dshift) { 3328fe8eb27SJed Brown if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);} 3338fe8eb27SJed Brown shell->dshift = shell->dshift_owned; 3348fe8eb27SJed Brown ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr); 3358fe8eb27SJed Brown } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);} 3368fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 3378fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 3388fe8eb27SJed Brown } else shell->vshift += a; 3398fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 340ef66eb69SBarry Smith PetscFunctionReturn(0); 341ef66eb69SBarry Smith } 342ef66eb69SBarry Smith 343ef66eb69SBarry Smith #undef __FUNCT__ 344ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell" 345f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 346ef66eb69SBarry Smith { 347ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3488fe8eb27SJed Brown PetscErrorCode ierr; 349b24ad042SBarry Smith 350ef66eb69SBarry Smith PetscFunctionBegin; 351f4df32b1SMatthew Knepley shell->vscale *= a; 3522205254eSKarl Rupp if (shell->dshift) { 3532205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 3542205254eSKarl Rupp } else shell->vshift *= a; 3558fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 3568fe8eb27SJed Brown PetscFunctionReturn(0); 35718be62a5SSatish Balay } 3588fe8eb27SJed Brown 3598fe8eb27SJed Brown #undef __FUNCT__ 3608fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell" 3618fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 3628fe8eb27SJed Brown { 3638fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 3648fe8eb27SJed Brown PetscErrorCode ierr; 3658fe8eb27SJed Brown 3668fe8eb27SJed Brown PetscFunctionBegin; 3678fe8eb27SJed Brown if (left) { 3688fe8eb27SJed Brown if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);} 3692205254eSKarl Rupp if (shell->left) { 3702205254eSKarl Rupp ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 3712205254eSKarl Rupp } else { 3728fe8eb27SJed Brown shell->left = shell->left_owned; 3738fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 37418be62a5SSatish Balay } 375ef66eb69SBarry Smith } 3768fe8eb27SJed Brown if (right) { 3778fe8eb27SJed Brown if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);} 3782205254eSKarl Rupp if (shell->right) { 3792205254eSKarl Rupp ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 3802205254eSKarl Rupp } else { 3818fe8eb27SJed Brown shell->right = shell->right_owned; 3828fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 3838fe8eb27SJed Brown } 3848fe8eb27SJed Brown } 3858fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 386ef66eb69SBarry Smith PetscFunctionReturn(0); 387ef66eb69SBarry Smith } 388ef66eb69SBarry Smith 389ef66eb69SBarry Smith #undef __FUNCT__ 390ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell" 391dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 392ef66eb69SBarry Smith { 393ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 394ef66eb69SBarry Smith 395ef66eb69SBarry Smith PetscFunctionBegin; 3968fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 397ef66eb69SBarry Smith shell->vshift = 0.0; 398ef66eb69SBarry Smith shell->vscale = 1.0; 3990298fd71SBarry Smith shell->dshift = NULL; 4000298fd71SBarry Smith shell->left = NULL; 4010298fd71SBarry Smith shell->right = NULL; 4027fb7d96aSJed Brown if (shell->mult) { 403ef66eb69SBarry Smith Y->ops->mult = shell->mult; 4040298fd71SBarry Smith shell->mult = NULL; 4057fb7d96aSJed Brown } 4067fb7d96aSJed Brown if (shell->multtranspose) { 40718be62a5SSatish Balay Y->ops->multtranspose = shell->multtranspose; 4080298fd71SBarry Smith shell->multtranspose = NULL; 4097fb7d96aSJed Brown } 4107fb7d96aSJed Brown if (shell->getdiagonal) { 41118be62a5SSatish Balay Y->ops->getdiagonal = shell->getdiagonal; 4120298fd71SBarry Smith shell->getdiagonal = NULL; 4137fb7d96aSJed Brown } 4147fb7d96aSJed Brown shell->usingscaled = PETSC_FALSE; 415ef66eb69SBarry Smith } 416ef66eb69SBarry Smith PetscFunctionReturn(0); 417ef66eb69SBarry Smith } 418ef66eb69SBarry Smith 41919fd82e9SBarry Smith extern PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*); 420b951964fSBarry Smith 42109dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 42220563c6bSBarry Smith 0, 42320563c6bSBarry Smith 0, 42420563c6bSBarry Smith 0, 42597304618SKris Buschelman /* 4*/ 0, 42620563c6bSBarry Smith 0, 427b951964fSBarry Smith 0, 428b951964fSBarry Smith 0, 429b951964fSBarry Smith 0, 430b951964fSBarry Smith 0, 43197304618SKris Buschelman /*10*/ 0, 432b951964fSBarry Smith 0, 433b951964fSBarry Smith 0, 434b951964fSBarry Smith 0, 435b951964fSBarry Smith 0, 43697304618SKris Buschelman /*15*/ 0, 437b951964fSBarry Smith 0, 438b951964fSBarry Smith 0, 4398fe8eb27SJed Brown MatDiagonalScale_Shell, 440b951964fSBarry Smith 0, 44197304618SKris Buschelman /*20*/ 0, 442ef66eb69SBarry Smith MatAssemblyEnd_Shell, 443b951964fSBarry Smith 0, 444b951964fSBarry Smith 0, 445d519adbfSMatthew Knepley /*24*/ 0, 446b951964fSBarry Smith 0, 447b951964fSBarry Smith 0, 448b951964fSBarry Smith 0, 449b951964fSBarry Smith 0, 450d519adbfSMatthew Knepley /*29*/ 0, 451b951964fSBarry Smith 0, 452273d9f13SBarry Smith 0, 453b951964fSBarry Smith 0, 454b951964fSBarry Smith 0, 455d519adbfSMatthew Knepley /*34*/ 0, 456b951964fSBarry Smith 0, 457b951964fSBarry Smith 0, 45809dc0095SBarry Smith 0, 45909dc0095SBarry Smith 0, 460d519adbfSMatthew Knepley /*39*/ 0, 46109dc0095SBarry Smith 0, 46209dc0095SBarry Smith 0, 46309dc0095SBarry Smith 0, 46409dc0095SBarry Smith 0, 465d519adbfSMatthew Knepley /*44*/ 0, 466ef66eb69SBarry Smith MatScale_Shell, 467ef66eb69SBarry Smith MatShift_Shell, 46809dc0095SBarry Smith 0, 46909dc0095SBarry Smith 0, 470f73d5cc4SBarry Smith /*49*/ 0, 47109dc0095SBarry Smith 0, 47209dc0095SBarry Smith 0, 47309dc0095SBarry Smith 0, 47409dc0095SBarry Smith 0, 475d519adbfSMatthew Knepley /*54*/ 0, 47609dc0095SBarry Smith 0, 47709dc0095SBarry Smith 0, 47809dc0095SBarry Smith 0, 47909dc0095SBarry Smith 0, 480d519adbfSMatthew Knepley /*59*/ 0, 481b9b97703SBarry Smith MatDestroy_Shell, 48209dc0095SBarry Smith 0, 483357abbc8SBarry Smith 0, 484273d9f13SBarry Smith 0, 485d519adbfSMatthew Knepley /*64*/ 0, 486273d9f13SBarry Smith 0, 487273d9f13SBarry Smith 0, 488273d9f13SBarry Smith 0, 489273d9f13SBarry Smith 0, 490d519adbfSMatthew Knepley /*69*/ 0, 491273d9f13SBarry Smith 0, 492c87e5d42SMatthew Knepley MatConvert_Shell, 493273d9f13SBarry Smith 0, 49497304618SKris Buschelman 0, 495d519adbfSMatthew Knepley /*74*/ 0, 49697304618SKris Buschelman 0, 49797304618SKris Buschelman 0, 49897304618SKris Buschelman 0, 49997304618SKris Buschelman 0, 500d519adbfSMatthew Knepley /*79*/ 0, 50197304618SKris Buschelman 0, 50297304618SKris Buschelman 0, 50397304618SKris Buschelman 0, 504865e5f61SKris Buschelman 0, 505d519adbfSMatthew Knepley /*84*/ 0, 506865e5f61SKris Buschelman 0, 507865e5f61SKris Buschelman 0, 508865e5f61SKris Buschelman 0, 509865e5f61SKris Buschelman 0, 510d519adbfSMatthew Knepley /*89*/ 0, 511865e5f61SKris Buschelman 0, 512865e5f61SKris Buschelman 0, 513865e5f61SKris Buschelman 0, 514865e5f61SKris Buschelman 0, 515d519adbfSMatthew Knepley /*94*/ 0, 516865e5f61SKris Buschelman 0, 517865e5f61SKris Buschelman 0, 5183964eb88SJed Brown 0, 5193964eb88SJed Brown 0, 5203964eb88SJed Brown /*99*/ 0, 5213964eb88SJed Brown 0, 5223964eb88SJed Brown 0, 5233964eb88SJed Brown 0, 5243964eb88SJed Brown 0, 5253964eb88SJed Brown /*104*/ 0, 5263964eb88SJed Brown 0, 5273964eb88SJed Brown 0, 5283964eb88SJed Brown 0, 5293964eb88SJed Brown 0, 5303964eb88SJed Brown /*109*/ 0, 5313964eb88SJed Brown 0, 5323964eb88SJed Brown 0, 5333964eb88SJed Brown 0, 5343964eb88SJed Brown 0, 5353964eb88SJed Brown /*114*/ 0, 5363964eb88SJed Brown 0, 5373964eb88SJed Brown 0, 5383964eb88SJed Brown 0, 5393964eb88SJed Brown 0, 5403964eb88SJed Brown /*119*/ 0, 5413964eb88SJed Brown 0, 5423964eb88SJed Brown 0, 5433964eb88SJed Brown 0, 5443964eb88SJed Brown 0, 5453964eb88SJed Brown /*124*/ 0, 5463964eb88SJed Brown 0, 5473964eb88SJed Brown 0, 5483964eb88SJed Brown 0, 5493964eb88SJed Brown 0, 5503964eb88SJed Brown /*129*/ 0, 5513964eb88SJed Brown 0, 5523964eb88SJed Brown 0, 5533964eb88SJed Brown 0, 5543964eb88SJed Brown 0, 5553964eb88SJed Brown /*134*/ 0, 5563964eb88SJed Brown 0, 5573964eb88SJed Brown 0, 5583964eb88SJed Brown 0, 5593964eb88SJed Brown 0, 5603964eb88SJed Brown /*139*/ 0, 561*f9426fe0SMark Adams 0, 5623964eb88SJed Brown 0 5633964eb88SJed Brown }; 564273d9f13SBarry Smith 5650bad9183SKris Buschelman /*MC 566fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 5670bad9183SKris Buschelman 5680bad9183SKris Buschelman Level: advanced 5690bad9183SKris Buschelman 5700bad9183SKris Buschelman .seealso: MatCreateShell 5710bad9183SKris Buschelman M*/ 5720bad9183SKris Buschelman 573711e205bSSatish Balay #undef __FUNCT__ 574711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell" 5758cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 576273d9f13SBarry Smith { 577273d9f13SBarry Smith Mat_Shell *b; 578dfbe8321SBarry Smith PetscErrorCode ierr; 579273d9f13SBarry Smith 580273d9f13SBarry Smith PetscFunctionBegin; 581273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 582273d9f13SBarry Smith 58338f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_Shell,&b);CHKERRQ(ierr); 584273d9f13SBarry Smith A->data = (void*)b; 585273d9f13SBarry Smith 58626283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 58726283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 588273d9f13SBarry Smith 589273d9f13SBarry Smith b->ctx = 0; 590ef66eb69SBarry Smith b->vshift = 0.0; 591ef66eb69SBarry Smith b->vscale = 1.0; 592ef66eb69SBarry Smith b->mult = 0; 59318be62a5SSatish Balay b->multtranspose = 0; 59418be62a5SSatish Balay b->getdiagonal = 0; 595273d9f13SBarry Smith A->assembled = PETSC_TRUE; 596210f0121SBarry Smith A->preallocated = PETSC_FALSE; 5972205254eSKarl Rupp 59817667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 599273d9f13SBarry Smith PetscFunctionReturn(0); 600273d9f13SBarry Smith } 601e51e0e81SBarry Smith 602711e205bSSatish Balay #undef __FUNCT__ 603711e205bSSatish Balay #define __FUNCT__ "MatCreateShell" 6044b828684SBarry Smith /*@C 605052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 606ff756334SLois Curfman McInnes private data storage format. 607e51e0e81SBarry Smith 608c7fcc2eaSBarry Smith Collective on MPI_Comm 609c7fcc2eaSBarry Smith 610e51e0e81SBarry Smith Input Parameters: 611c7fcc2eaSBarry Smith + comm - MPI communicator 612c7fcc2eaSBarry Smith . m - number of local rows (must be given) 613c7fcc2eaSBarry Smith . n - number of local columns (must be given) 614c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 615c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 616c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 617e51e0e81SBarry Smith 618ff756334SLois Curfman McInnes Output Parameter: 61944cd7ae7SLois Curfman McInnes . A - the matrix 620e51e0e81SBarry Smith 621ff2fd236SBarry Smith Level: advanced 622ff2fd236SBarry Smith 623f39d1f56SLois Curfman McInnes Usage: 6247b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 625f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 626c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 627f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 628f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 629f39d1f56SLois Curfman McInnes 630ff756334SLois Curfman McInnes Notes: 631ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 632ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 633ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 634e51e0e81SBarry Smith 635f38a8d56SBarry Smith Fortran Notes: The context can only be an integer or a PetscObject 636f38a8d56SBarry Smith unfortunately it cannot be a Fortran array or derived type. 637f38a8d56SBarry Smith 638f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 639f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 640645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 641645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 642645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 643645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 644645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 645645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 646645985a0SLois Curfman McInnes For example, 647f39d1f56SLois Curfman McInnes 648f39d1f56SLois Curfman McInnes $ 649f39d1f56SLois Curfman McInnes $ Vec x, y 6507b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 651645985a0SLois Curfman McInnes $ Mat A 652f39d1f56SLois Curfman McInnes $ 653c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 654c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 655f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 656c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 657c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 658c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 659645985a0SLois Curfman McInnes $ MatMult(A,x,y); 660645985a0SLois Curfman McInnes $ MatDestroy(A); 661f39d1f56SLois Curfman McInnes $ VecDestroy(y); VecDestroy(x); 662645985a0SLois Curfman McInnes $ 663e51e0e81SBarry Smith 6640b627109SLois Curfman McInnes .keywords: matrix, shell, create 6650b627109SLois Curfman McInnes 666ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext() 667e51e0e81SBarry Smith @*/ 6687087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 669e51e0e81SBarry Smith { 670dfbe8321SBarry Smith PetscErrorCode ierr; 671ed3cc1f0SBarry Smith 6723a40ed3dSBarry Smith PetscFunctionBegin; 673f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 674f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 675273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 676273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 6770fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 678273d9f13SBarry Smith PetscFunctionReturn(0); 679c7fcc2eaSBarry Smith } 680c7fcc2eaSBarry Smith 681711e205bSSatish Balay #undef __FUNCT__ 682711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext" 683c6866cfdSSatish Balay /*@ 684273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 685c7fcc2eaSBarry Smith 6863f9fe445SBarry Smith Logically Collective on Mat 687c7fcc2eaSBarry Smith 688273d9f13SBarry Smith Input Parameters: 689273d9f13SBarry Smith + mat - the shell matrix 690273d9f13SBarry Smith - ctx - the context 691273d9f13SBarry Smith 692273d9f13SBarry Smith Level: advanced 693273d9f13SBarry Smith 694f38a8d56SBarry Smith Fortran Notes: The context can only be an integer or a PetscObject 695f38a8d56SBarry Smith unfortunately it cannot be a Fortran array or derived type. 696273d9f13SBarry Smith 697273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 6980bc0a719SBarry Smith @*/ 6997087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 700273d9f13SBarry Smith { 701273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 702dfbe8321SBarry Smith PetscErrorCode ierr; 703ace3abfcSBarry Smith PetscBool flg; 704273d9f13SBarry Smith 705273d9f13SBarry Smith PetscFunctionBegin; 7060700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 707251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 708273d9f13SBarry Smith if (flg) { 709273d9f13SBarry Smith shell->ctx = ctx; 710ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix"); 7113a40ed3dSBarry Smith PetscFunctionReturn(0); 712e51e0e81SBarry Smith } 713e51e0e81SBarry Smith 714711e205bSSatish Balay #undef __FUNCT__ 715711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation" 716c16cb8f2SBarry Smith /*@C 7173a3eedf2SBarry Smith MatShellSetOperation - Allows user to set a matrix operation for 7183a3eedf2SBarry Smith a shell matrix. 719e51e0e81SBarry Smith 7203f9fe445SBarry Smith Logically Collective on Mat 721fee21e36SBarry Smith 722c7fcc2eaSBarry Smith Input Parameters: 723c7fcc2eaSBarry Smith + mat - the shell matrix 724c7fcc2eaSBarry Smith . op - the name of the operation 725c7fcc2eaSBarry Smith - f - the function that provides the operation. 726c7fcc2eaSBarry Smith 72715091d37SBarry Smith Level: advanced 72815091d37SBarry Smith 729fae171e0SBarry Smith Usage: 730e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 731f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 732c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 7330b627109SLois Curfman McInnes 734a62d957aSLois Curfman McInnes Notes: 735e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 7361c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 737a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 7381c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 739a62d957aSLois Curfman McInnes 74025296bd5SBarry Smith All user-provided functions (execept for MATOP_DESTROY) should have the same calling 741deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 742deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 743deebb3c3SLois Curfman McInnes routines, e.g., 744a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 745a62d957aSLois Curfman McInnes 7464aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 7474aa34b0aSBarry Smith nonzero on failure. 7484aa34b0aSBarry Smith 749a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 750a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 751a62d957aSLois Curfman McInnes set by MatCreateShell(). 752a62d957aSLois Curfman McInnes 753501d9185SBarry Smith Fortran Notes: For MatGetVecs() the user code should check if the input left or right matrix is -1 and in that case not 754501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 755501d9185SBarry Smith 756a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation 757a62d957aSLois Curfman McInnes 758ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext() 759e51e0e81SBarry Smith @*/ 7607087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 761e51e0e81SBarry Smith { 762dfbe8321SBarry Smith PetscErrorCode ierr; 763ace3abfcSBarry Smith PetscBool flg; 764273d9f13SBarry Smith 7653a40ed3dSBarry Smith PetscFunctionBegin; 7660700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 7675edf6516SJed Brown switch (op) { 7685edf6516SJed Brown case MATOP_DESTROY: 769251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 770273d9f13SBarry Smith if (flg) { 771a62d957aSLois Curfman McInnes Mat_Shell *shell = (Mat_Shell*)mat->data; 7726849ba73SBarry Smith shell->destroy = (PetscErrorCode (*)(Mat))f; 7736849ba73SBarry Smith } else mat->ops->destroy = (PetscErrorCode (*)(Mat))f; 7745edf6516SJed Brown break; 7755edf6516SJed Brown case MATOP_VIEW: 7765edf6516SJed Brown mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 7775edf6516SJed Brown break; 7785edf6516SJed Brown case MATOP_MULT: 7795edf6516SJed Brown mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 7805edf6516SJed Brown if (!mat->ops->multadd) mat->ops->multadd = MatMultAdd_Shell; 7815edf6516SJed Brown break; 7825edf6516SJed Brown case MATOP_MULT_TRANSPOSE: 7835edf6516SJed Brown mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 7845edf6516SJed Brown if (!mat->ops->multtransposeadd) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell; 7855edf6516SJed Brown break; 7865edf6516SJed Brown default: 7875edf6516SJed Brown (((void(**)(void))mat->ops)[op]) = f; 788a62d957aSLois Curfman McInnes } 7893a40ed3dSBarry Smith PetscFunctionReturn(0); 790e51e0e81SBarry Smith } 791f0479e8cSBarry Smith 792711e205bSSatish Balay #undef __FUNCT__ 793711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation" 794d4bb536fSBarry Smith /*@C 795d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 796d4bb536fSBarry Smith 797c7fcc2eaSBarry Smith Not Collective 798c7fcc2eaSBarry Smith 799d4bb536fSBarry Smith Input Parameters: 800c7fcc2eaSBarry Smith + mat - the shell matrix 801c7fcc2eaSBarry Smith - op - the name of the operation 802d4bb536fSBarry Smith 803d4bb536fSBarry Smith Output Parameter: 804d4bb536fSBarry Smith . f - the function that provides the operation. 805d4bb536fSBarry Smith 80615091d37SBarry Smith Level: advanced 80715091d37SBarry Smith 808d4bb536fSBarry Smith Notes: 809e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 810d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 811d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 812d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 813d4bb536fSBarry Smith 814d4bb536fSBarry Smith All user-provided functions have the same calling 815d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 816d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 817d4bb536fSBarry Smith routines, e.g., 818d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 819d4bb536fSBarry Smith 820d4bb536fSBarry Smith Within each user-defined routine, the user should call 821d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 822d4bb536fSBarry Smith set by MatCreateShell(). 823d4bb536fSBarry Smith 824d4bb536fSBarry Smith .keywords: matrix, shell, set, operation 825d4bb536fSBarry Smith 826ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 827d4bb536fSBarry Smith @*/ 8287087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 829d4bb536fSBarry Smith { 830dfbe8321SBarry Smith PetscErrorCode ierr; 831ace3abfcSBarry Smith PetscBool flg; 832273d9f13SBarry Smith 8333a40ed3dSBarry Smith PetscFunctionBegin; 8340700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 835d4bb536fSBarry Smith if (op == MATOP_DESTROY) { 836251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 837273d9f13SBarry Smith if (flg) { 838d4bb536fSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 839c134de8dSSatish Balay *f = (void (*)(void))shell->destroy; 840c7fcc2eaSBarry Smith } else { 841c134de8dSSatish Balay *f = (void (*)(void))mat->ops->destroy; 842d4bb536fSBarry Smith } 843c7fcc2eaSBarry Smith } else if (op == MATOP_VIEW) { 844c134de8dSSatish Balay *f = (void (*)(void))mat->ops->view; 845c7fcc2eaSBarry Smith } else { 846c134de8dSSatish Balay *f = (((void (**)(void))mat->ops)[op]); 847d4bb536fSBarry Smith } 8483a40ed3dSBarry Smith PetscFunctionReturn(0); 849d4bb536fSBarry Smith } 850d4bb536fSBarry Smith 851