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*/ 9af0996ceSBarry Smith #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 188daf670e6SBarry Smith Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this 189daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 190a62d957aSLois Curfman McInnes 191a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context 192a62d957aSLois Curfman McInnes 193ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 1940bc0a719SBarry Smith @*/ 1958fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 196b4fd4287SBarry Smith { 197dfbe8321SBarry Smith PetscErrorCode ierr; 198ace3abfcSBarry Smith PetscBool flg; 199273d9f13SBarry Smith 2003a40ed3dSBarry Smith PetscFunctionBegin; 2010700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2024482741eSBarry Smith PetscValidPointer(ctx,2); 203251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 204940183f3SBarry Smith if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 205ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix"); 2063a40ed3dSBarry Smith PetscFunctionReturn(0); 207b4fd4287SBarry Smith } 208b4fd4287SBarry Smith 209711e205bSSatish Balay #undef __FUNCT__ 210711e205bSSatish Balay #define __FUNCT__ "MatDestroy_Shell" 211dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 212e51e0e81SBarry Smith { 213dfbe8321SBarry Smith PetscErrorCode ierr; 214bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 215ed3cc1f0SBarry Smith 2163a40ed3dSBarry Smith PetscFunctionBegin; 217bf0cc555SLisandro Dalcin if (shell->destroy) { 218bf0cc555SLisandro Dalcin ierr = (*shell->destroy)(mat);CHKERRQ(ierr); 219bf0cc555SLisandro Dalcin } 2208fe8eb27SJed Brown ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr); 2218fe8eb27SJed Brown ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr); 2228fe8eb27SJed Brown ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr); 2238fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 2248fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 2255edf6516SJed Brown ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr); 2265edf6516SJed Brown ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr); 227bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 2283a40ed3dSBarry Smith PetscFunctionReturn(0); 229e51e0e81SBarry Smith } 230e51e0e81SBarry Smith 231711e205bSSatish Balay #undef __FUNCT__ 232ef66eb69SBarry Smith #define __FUNCT__ "MatMult_Shell" 233dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 234ef66eb69SBarry Smith { 235ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 236dfbe8321SBarry Smith PetscErrorCode ierr; 23725578ef6SJed Brown Vec xx; 238ef66eb69SBarry Smith 239ef66eb69SBarry Smith PetscFunctionBegin; 2408fe8eb27SJed Brown ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr); 2418fe8eb27SJed Brown ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr); 2428fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2438fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 244ef66eb69SBarry Smith PetscFunctionReturn(0); 245ef66eb69SBarry Smith } 246ef66eb69SBarry Smith 247ef66eb69SBarry Smith #undef __FUNCT__ 2485edf6516SJed Brown #define __FUNCT__ "MatMultAdd_Shell" 2495edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 2505edf6516SJed Brown { 2515edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2525edf6516SJed Brown PetscErrorCode ierr; 2535edf6516SJed Brown 2545edf6516SJed Brown PetscFunctionBegin; 2555edf6516SJed Brown if (y == z) { 2565edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 2575edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 258b8430d49SBarry Smith ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr); 2595edf6516SJed Brown } else { 2605edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 2615edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 2625edf6516SJed Brown } 2635edf6516SJed Brown PetscFunctionReturn(0); 2645edf6516SJed Brown } 2655edf6516SJed Brown 2665edf6516SJed Brown #undef __FUNCT__ 26718be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell" 26818be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 26918be62a5SSatish Balay { 27018be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 27118be62a5SSatish Balay PetscErrorCode ierr; 27225578ef6SJed Brown Vec xx; 27318be62a5SSatish Balay 27418be62a5SSatish Balay PetscFunctionBegin; 2758fe8eb27SJed Brown ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr); 2768fe8eb27SJed Brown ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr); 2778fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2788fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 27918be62a5SSatish Balay PetscFunctionReturn(0); 28018be62a5SSatish Balay } 28118be62a5SSatish Balay 28218be62a5SSatish Balay #undef __FUNCT__ 2835edf6516SJed Brown #define __FUNCT__ "MatMultTransposeAdd_Shell" 2845edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 2855edf6516SJed Brown { 2865edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2875edf6516SJed Brown PetscErrorCode ierr; 2885edf6516SJed Brown 2895edf6516SJed Brown PetscFunctionBegin; 2905edf6516SJed Brown if (y == z) { 2915edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 2925edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 2935edf6516SJed Brown ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr); 2945edf6516SJed Brown } else { 2955edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 2965edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 2975edf6516SJed Brown } 2985edf6516SJed Brown PetscFunctionReturn(0); 2995edf6516SJed Brown } 3005edf6516SJed Brown 3015edf6516SJed Brown #undef __FUNCT__ 30218be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell" 30318be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 30418be62a5SSatish Balay { 30518be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 30618be62a5SSatish Balay PetscErrorCode ierr; 30718be62a5SSatish Balay 30818be62a5SSatish Balay PetscFunctionBegin; 30918be62a5SSatish Balay ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr); 31018be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 3118fe8eb27SJed Brown if (shell->dshift) { 3128fe8eb27SJed Brown ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr); 3138fe8eb27SJed Brown } else { 31418be62a5SSatish Balay ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 31518be62a5SSatish Balay } 3168fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 3178fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 31818be62a5SSatish Balay PetscFunctionReturn(0); 31918be62a5SSatish Balay } 32018be62a5SSatish Balay 32118be62a5SSatish Balay #undef __FUNCT__ 322ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell" 323f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 324ef66eb69SBarry Smith { 325ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3268fe8eb27SJed Brown PetscErrorCode ierr; 327b24ad042SBarry Smith 328ef66eb69SBarry Smith PetscFunctionBegin; 3298fe8eb27SJed Brown if (shell->left || shell->right || shell->dshift) { 3308fe8eb27SJed Brown if (!shell->dshift) { 3318fe8eb27SJed Brown if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);} 3328fe8eb27SJed Brown shell->dshift = shell->dshift_owned; 3338fe8eb27SJed Brown ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr); 3348fe8eb27SJed Brown } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);} 3358fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 3368fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 3378fe8eb27SJed Brown } else shell->vshift += a; 3388fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 339ef66eb69SBarry Smith PetscFunctionReturn(0); 340ef66eb69SBarry Smith } 341ef66eb69SBarry Smith 342ef66eb69SBarry Smith #undef __FUNCT__ 343ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell" 344f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 345ef66eb69SBarry Smith { 346ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3478fe8eb27SJed Brown PetscErrorCode ierr; 348b24ad042SBarry Smith 349ef66eb69SBarry Smith PetscFunctionBegin; 350f4df32b1SMatthew Knepley shell->vscale *= a; 3512205254eSKarl Rupp if (shell->dshift) { 3522205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 3532205254eSKarl Rupp } else shell->vshift *= a; 3548fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 3558fe8eb27SJed Brown PetscFunctionReturn(0); 35618be62a5SSatish Balay } 3578fe8eb27SJed Brown 3588fe8eb27SJed Brown #undef __FUNCT__ 3598fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell" 3608fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 3618fe8eb27SJed Brown { 3628fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 3638fe8eb27SJed Brown PetscErrorCode ierr; 3648fe8eb27SJed Brown 3658fe8eb27SJed Brown PetscFunctionBegin; 3668fe8eb27SJed Brown if (left) { 3678fe8eb27SJed Brown if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);} 3682205254eSKarl Rupp if (shell->left) { 3692205254eSKarl Rupp ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 3702205254eSKarl Rupp } else { 3718fe8eb27SJed Brown shell->left = shell->left_owned; 3728fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 37318be62a5SSatish Balay } 374ef66eb69SBarry Smith } 3758fe8eb27SJed Brown if (right) { 3768fe8eb27SJed Brown if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);} 3772205254eSKarl Rupp if (shell->right) { 3782205254eSKarl Rupp ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 3792205254eSKarl Rupp } else { 3808fe8eb27SJed Brown shell->right = shell->right_owned; 3818fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 3828fe8eb27SJed Brown } 3838fe8eb27SJed Brown } 3848fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 385ef66eb69SBarry Smith PetscFunctionReturn(0); 386ef66eb69SBarry Smith } 387ef66eb69SBarry Smith 388ef66eb69SBarry Smith #undef __FUNCT__ 389ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell" 390dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 391ef66eb69SBarry Smith { 392ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 393ef66eb69SBarry Smith 394ef66eb69SBarry Smith PetscFunctionBegin; 3958fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 396ef66eb69SBarry Smith shell->vshift = 0.0; 397ef66eb69SBarry Smith shell->vscale = 1.0; 3980298fd71SBarry Smith shell->dshift = NULL; 3990298fd71SBarry Smith shell->left = NULL; 4000298fd71SBarry Smith shell->right = NULL; 4017fb7d96aSJed Brown if (shell->mult) { 402ef66eb69SBarry Smith Y->ops->mult = shell->mult; 4030298fd71SBarry Smith shell->mult = NULL; 4047fb7d96aSJed Brown } 4057fb7d96aSJed Brown if (shell->multtranspose) { 40618be62a5SSatish Balay Y->ops->multtranspose = shell->multtranspose; 4070298fd71SBarry Smith shell->multtranspose = NULL; 4087fb7d96aSJed Brown } 4097fb7d96aSJed Brown if (shell->getdiagonal) { 41018be62a5SSatish Balay Y->ops->getdiagonal = shell->getdiagonal; 4110298fd71SBarry Smith shell->getdiagonal = NULL; 4127fb7d96aSJed Brown } 4137fb7d96aSJed Brown shell->usingscaled = PETSC_FALSE; 414ef66eb69SBarry Smith } 415ef66eb69SBarry Smith PetscFunctionReturn(0); 416ef66eb69SBarry Smith } 417ef66eb69SBarry Smith 41819fd82e9SBarry Smith extern PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*); 419b951964fSBarry Smith 420*3b49f96aSBarry Smith #undef __FUNCT__ 421*3b49f96aSBarry Smith #define __FUNCT__ "MatMissingDiagonal_Shell" 422*3b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 423*3b49f96aSBarry Smith { 424*3b49f96aSBarry Smith PetscFunctionBegin; 425*3b49f96aSBarry Smith *missing = PETSC_FALSE; 426*3b49f96aSBarry Smith PetscFunctionReturn(0); 427*3b49f96aSBarry Smith } 428*3b49f96aSBarry Smith 42909dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 43020563c6bSBarry Smith 0, 43120563c6bSBarry Smith 0, 43220563c6bSBarry Smith 0, 43397304618SKris Buschelman /* 4*/ 0, 43420563c6bSBarry Smith 0, 435b951964fSBarry Smith 0, 436b951964fSBarry Smith 0, 437b951964fSBarry Smith 0, 438b951964fSBarry Smith 0, 43997304618SKris Buschelman /*10*/ 0, 440b951964fSBarry Smith 0, 441b951964fSBarry Smith 0, 442b951964fSBarry Smith 0, 443b951964fSBarry Smith 0, 44497304618SKris Buschelman /*15*/ 0, 445b951964fSBarry Smith 0, 446b951964fSBarry Smith 0, 4478fe8eb27SJed Brown MatDiagonalScale_Shell, 448b951964fSBarry Smith 0, 44997304618SKris Buschelman /*20*/ 0, 450ef66eb69SBarry Smith MatAssemblyEnd_Shell, 451b951964fSBarry Smith 0, 452b951964fSBarry Smith 0, 453d519adbfSMatthew Knepley /*24*/ 0, 454b951964fSBarry Smith 0, 455b951964fSBarry Smith 0, 456b951964fSBarry Smith 0, 457b951964fSBarry Smith 0, 458d519adbfSMatthew Knepley /*29*/ 0, 459b951964fSBarry Smith 0, 460273d9f13SBarry Smith 0, 461b951964fSBarry Smith 0, 462b951964fSBarry Smith 0, 463d519adbfSMatthew Knepley /*34*/ 0, 464b951964fSBarry Smith 0, 465b951964fSBarry Smith 0, 46609dc0095SBarry Smith 0, 46709dc0095SBarry Smith 0, 468d519adbfSMatthew Knepley /*39*/ 0, 46909dc0095SBarry Smith 0, 47009dc0095SBarry Smith 0, 47109dc0095SBarry Smith 0, 47209dc0095SBarry Smith 0, 473d519adbfSMatthew Knepley /*44*/ 0, 474ef66eb69SBarry Smith MatScale_Shell, 475ef66eb69SBarry Smith MatShift_Shell, 47609dc0095SBarry Smith 0, 47709dc0095SBarry Smith 0, 478f73d5cc4SBarry Smith /*49*/ 0, 47909dc0095SBarry Smith 0, 48009dc0095SBarry Smith 0, 48109dc0095SBarry Smith 0, 48209dc0095SBarry Smith 0, 483d519adbfSMatthew Knepley /*54*/ 0, 48409dc0095SBarry Smith 0, 48509dc0095SBarry Smith 0, 48609dc0095SBarry Smith 0, 48709dc0095SBarry Smith 0, 488d519adbfSMatthew Knepley /*59*/ 0, 489b9b97703SBarry Smith MatDestroy_Shell, 49009dc0095SBarry Smith 0, 491357abbc8SBarry Smith 0, 492273d9f13SBarry Smith 0, 493d519adbfSMatthew Knepley /*64*/ 0, 494273d9f13SBarry Smith 0, 495273d9f13SBarry Smith 0, 496273d9f13SBarry Smith 0, 497273d9f13SBarry Smith 0, 498d519adbfSMatthew Knepley /*69*/ 0, 499273d9f13SBarry Smith 0, 500c87e5d42SMatthew Knepley MatConvert_Shell, 501273d9f13SBarry Smith 0, 50297304618SKris Buschelman 0, 503d519adbfSMatthew Knepley /*74*/ 0, 50497304618SKris Buschelman 0, 50597304618SKris Buschelman 0, 50697304618SKris Buschelman 0, 50797304618SKris Buschelman 0, 508d519adbfSMatthew Knepley /*79*/ 0, 50997304618SKris Buschelman 0, 51097304618SKris Buschelman 0, 51197304618SKris Buschelman 0, 512865e5f61SKris Buschelman 0, 513d519adbfSMatthew Knepley /*84*/ 0, 514865e5f61SKris Buschelman 0, 515865e5f61SKris Buschelman 0, 516865e5f61SKris Buschelman 0, 517865e5f61SKris Buschelman 0, 518d519adbfSMatthew Knepley /*89*/ 0, 519865e5f61SKris Buschelman 0, 520865e5f61SKris Buschelman 0, 521865e5f61SKris Buschelman 0, 522865e5f61SKris Buschelman 0, 523d519adbfSMatthew Knepley /*94*/ 0, 524865e5f61SKris Buschelman 0, 525865e5f61SKris Buschelman 0, 5263964eb88SJed Brown 0, 5273964eb88SJed Brown 0, 5283964eb88SJed Brown /*99*/ 0, 5293964eb88SJed Brown 0, 5303964eb88SJed Brown 0, 5313964eb88SJed Brown 0, 5323964eb88SJed Brown 0, 5333964eb88SJed Brown /*104*/ 0, 5343964eb88SJed Brown 0, 5353964eb88SJed Brown 0, 5363964eb88SJed Brown 0, 5373964eb88SJed Brown 0, 5383964eb88SJed Brown /*109*/ 0, 5393964eb88SJed Brown 0, 5403964eb88SJed Brown 0, 5413964eb88SJed Brown 0, 542*3b49f96aSBarry Smith MatMissingDiagonal_Shell, 5433964eb88SJed Brown /*114*/ 0, 5443964eb88SJed Brown 0, 5453964eb88SJed Brown 0, 5463964eb88SJed Brown 0, 5473964eb88SJed Brown 0, 5483964eb88SJed Brown /*119*/ 0, 5493964eb88SJed Brown 0, 5503964eb88SJed Brown 0, 5513964eb88SJed Brown 0, 5523964eb88SJed Brown 0, 5533964eb88SJed Brown /*124*/ 0, 5543964eb88SJed Brown 0, 5553964eb88SJed Brown 0, 5563964eb88SJed Brown 0, 5573964eb88SJed Brown 0, 5583964eb88SJed Brown /*129*/ 0, 5593964eb88SJed Brown 0, 5603964eb88SJed Brown 0, 5613964eb88SJed Brown 0, 5623964eb88SJed Brown 0, 5633964eb88SJed Brown /*134*/ 0, 5643964eb88SJed Brown 0, 5653964eb88SJed Brown 0, 5663964eb88SJed Brown 0, 5673964eb88SJed Brown 0, 5683964eb88SJed Brown /*139*/ 0, 569f9426fe0SMark Adams 0, 5703964eb88SJed Brown 0 5713964eb88SJed Brown }; 572273d9f13SBarry Smith 5730bad9183SKris Buschelman /*MC 574fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 5750bad9183SKris Buschelman 5760bad9183SKris Buschelman Level: advanced 5770bad9183SKris Buschelman 5780bad9183SKris Buschelman .seealso: MatCreateShell 5790bad9183SKris Buschelman M*/ 5800bad9183SKris Buschelman 581711e205bSSatish Balay #undef __FUNCT__ 582711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell" 5838cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 584273d9f13SBarry Smith { 585273d9f13SBarry Smith Mat_Shell *b; 586dfbe8321SBarry Smith PetscErrorCode ierr; 587273d9f13SBarry Smith 588273d9f13SBarry Smith PetscFunctionBegin; 589273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 590273d9f13SBarry Smith 591b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 592273d9f13SBarry Smith A->data = (void*)b; 593273d9f13SBarry Smith 59426283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 59526283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 596273d9f13SBarry Smith 597273d9f13SBarry Smith b->ctx = 0; 598ef66eb69SBarry Smith b->vshift = 0.0; 599ef66eb69SBarry Smith b->vscale = 1.0; 600ef66eb69SBarry Smith b->mult = 0; 60118be62a5SSatish Balay b->multtranspose = 0; 60218be62a5SSatish Balay b->getdiagonal = 0; 603273d9f13SBarry Smith A->assembled = PETSC_TRUE; 604210f0121SBarry Smith A->preallocated = PETSC_FALSE; 6052205254eSKarl Rupp 60617667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 607273d9f13SBarry Smith PetscFunctionReturn(0); 608273d9f13SBarry Smith } 609e51e0e81SBarry Smith 610711e205bSSatish Balay #undef __FUNCT__ 611711e205bSSatish Balay #define __FUNCT__ "MatCreateShell" 6124b828684SBarry Smith /*@C 613052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 614ff756334SLois Curfman McInnes private data storage format. 615e51e0e81SBarry Smith 616c7fcc2eaSBarry Smith Collective on MPI_Comm 617c7fcc2eaSBarry Smith 618e51e0e81SBarry Smith Input Parameters: 619c7fcc2eaSBarry Smith + comm - MPI communicator 620c7fcc2eaSBarry Smith . m - number of local rows (must be given) 621c7fcc2eaSBarry Smith . n - number of local columns (must be given) 622c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 623c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 624c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 625e51e0e81SBarry Smith 626ff756334SLois Curfman McInnes Output Parameter: 62744cd7ae7SLois Curfman McInnes . A - the matrix 628e51e0e81SBarry Smith 629ff2fd236SBarry Smith Level: advanced 630ff2fd236SBarry Smith 631f39d1f56SLois Curfman McInnes Usage: 6327b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 633f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 634c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 635f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 636f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 637f39d1f56SLois Curfman McInnes 638ff756334SLois Curfman McInnes Notes: 639ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 640ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 641ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 642e51e0e81SBarry Smith 643daf670e6SBarry Smith Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this 644daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 645daf670e6SBarry Smith in as the ctx argument. 646f38a8d56SBarry Smith 647f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 648f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 649645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 650645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 651645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 652645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 653645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 654645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 655645985a0SLois Curfman McInnes For example, 656f39d1f56SLois Curfman McInnes 657f39d1f56SLois Curfman McInnes $ 658f39d1f56SLois Curfman McInnes $ Vec x, y 6597b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 660645985a0SLois Curfman McInnes $ Mat A 661f39d1f56SLois Curfman McInnes $ 662c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 663c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 664f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 665c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 666c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 667c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 668645985a0SLois Curfman McInnes $ MatMult(A,x,y); 669645985a0SLois Curfman McInnes $ MatDestroy(A); 670f39d1f56SLois Curfman McInnes $ VecDestroy(y); VecDestroy(x); 671645985a0SLois Curfman McInnes $ 672e51e0e81SBarry Smith 6730b627109SLois Curfman McInnes .keywords: matrix, shell, create 6740b627109SLois Curfman McInnes 675ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext() 676e51e0e81SBarry Smith @*/ 6777087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 678e51e0e81SBarry Smith { 679dfbe8321SBarry Smith PetscErrorCode ierr; 680ed3cc1f0SBarry Smith 6813a40ed3dSBarry Smith PetscFunctionBegin; 682f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 683f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 684273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 685273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 6860fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 687273d9f13SBarry Smith PetscFunctionReturn(0); 688c7fcc2eaSBarry Smith } 689c7fcc2eaSBarry Smith 690711e205bSSatish Balay #undef __FUNCT__ 691711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext" 692c6866cfdSSatish Balay /*@ 693273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 694c7fcc2eaSBarry Smith 6953f9fe445SBarry Smith Logically Collective on Mat 696c7fcc2eaSBarry Smith 697273d9f13SBarry Smith Input Parameters: 698273d9f13SBarry Smith + mat - the shell matrix 699273d9f13SBarry Smith - ctx - the context 700273d9f13SBarry Smith 701273d9f13SBarry Smith Level: advanced 702273d9f13SBarry Smith 703daf670e6SBarry Smith Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this 704daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 705273d9f13SBarry Smith 706273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 7070bc0a719SBarry Smith @*/ 7087087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 709273d9f13SBarry Smith { 710273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 711dfbe8321SBarry Smith PetscErrorCode ierr; 712ace3abfcSBarry Smith PetscBool flg; 713273d9f13SBarry Smith 714273d9f13SBarry Smith PetscFunctionBegin; 7150700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 716251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 717273d9f13SBarry Smith if (flg) { 718273d9f13SBarry Smith shell->ctx = ctx; 719ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix"); 7203a40ed3dSBarry Smith PetscFunctionReturn(0); 721e51e0e81SBarry Smith } 722e51e0e81SBarry Smith 723711e205bSSatish Balay #undef __FUNCT__ 724711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation" 725c16cb8f2SBarry Smith /*@C 7263a3eedf2SBarry Smith MatShellSetOperation - Allows user to set a matrix operation for 7273a3eedf2SBarry Smith a shell matrix. 728e51e0e81SBarry Smith 7293f9fe445SBarry Smith Logically Collective on Mat 730fee21e36SBarry Smith 731c7fcc2eaSBarry Smith Input Parameters: 732c7fcc2eaSBarry Smith + mat - the shell matrix 733c7fcc2eaSBarry Smith . op - the name of the operation 734c7fcc2eaSBarry Smith - f - the function that provides the operation. 735c7fcc2eaSBarry Smith 73615091d37SBarry Smith Level: advanced 73715091d37SBarry Smith 738fae171e0SBarry Smith Usage: 739e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 740f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 741c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 7420b627109SLois Curfman McInnes 743a62d957aSLois Curfman McInnes Notes: 744e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 7451c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 746a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 7471c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 748a62d957aSLois Curfman McInnes 74925296bd5SBarry Smith All user-provided functions (execept for MATOP_DESTROY) should have the same calling 750deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 751deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 752deebb3c3SLois Curfman McInnes routines, e.g., 753a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 754a62d957aSLois Curfman McInnes 7554aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 7564aa34b0aSBarry Smith nonzero on failure. 7574aa34b0aSBarry Smith 758a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 759a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 760a62d957aSLois Curfman McInnes set by MatCreateShell(). 761a62d957aSLois Curfman McInnes 7622a7a6963SBarry Smith Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 763501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 764501d9185SBarry Smith 765a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation 766a62d957aSLois Curfman McInnes 767ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext() 768e51e0e81SBarry Smith @*/ 7697087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 770e51e0e81SBarry Smith { 771dfbe8321SBarry Smith PetscErrorCode ierr; 772ace3abfcSBarry Smith PetscBool flg; 773273d9f13SBarry Smith 7743a40ed3dSBarry Smith PetscFunctionBegin; 7750700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 7765edf6516SJed Brown switch (op) { 7775edf6516SJed Brown case MATOP_DESTROY: 778251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 779273d9f13SBarry Smith if (flg) { 780a62d957aSLois Curfman McInnes Mat_Shell *shell = (Mat_Shell*)mat->data; 7816849ba73SBarry Smith shell->destroy = (PetscErrorCode (*)(Mat))f; 7826849ba73SBarry Smith } else mat->ops->destroy = (PetscErrorCode (*)(Mat))f; 7835edf6516SJed Brown break; 7845edf6516SJed Brown case MATOP_VIEW: 7855edf6516SJed Brown mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 7865edf6516SJed Brown break; 7875edf6516SJed Brown case MATOP_MULT: 7885edf6516SJed Brown mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 7895edf6516SJed Brown if (!mat->ops->multadd) mat->ops->multadd = MatMultAdd_Shell; 7905edf6516SJed Brown break; 7915edf6516SJed Brown case MATOP_MULT_TRANSPOSE: 7925edf6516SJed Brown mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 7935edf6516SJed Brown if (!mat->ops->multtransposeadd) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell; 7945edf6516SJed Brown break; 7955edf6516SJed Brown default: 7965edf6516SJed Brown (((void(**)(void))mat->ops)[op]) = f; 797a62d957aSLois Curfman McInnes } 7983a40ed3dSBarry Smith PetscFunctionReturn(0); 799e51e0e81SBarry Smith } 800f0479e8cSBarry Smith 801711e205bSSatish Balay #undef __FUNCT__ 802711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation" 803d4bb536fSBarry Smith /*@C 804d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 805d4bb536fSBarry Smith 806c7fcc2eaSBarry Smith Not Collective 807c7fcc2eaSBarry Smith 808d4bb536fSBarry Smith Input Parameters: 809c7fcc2eaSBarry Smith + mat - the shell matrix 810c7fcc2eaSBarry Smith - op - the name of the operation 811d4bb536fSBarry Smith 812d4bb536fSBarry Smith Output Parameter: 813d4bb536fSBarry Smith . f - the function that provides the operation. 814d4bb536fSBarry Smith 81515091d37SBarry Smith Level: advanced 81615091d37SBarry Smith 817d4bb536fSBarry Smith Notes: 818e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 819d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 820d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 821d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 822d4bb536fSBarry Smith 823d4bb536fSBarry Smith All user-provided functions have the same calling 824d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 825d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 826d4bb536fSBarry Smith routines, e.g., 827d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 828d4bb536fSBarry Smith 829d4bb536fSBarry Smith Within each user-defined routine, the user should call 830d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 831d4bb536fSBarry Smith set by MatCreateShell(). 832d4bb536fSBarry Smith 833d4bb536fSBarry Smith .keywords: matrix, shell, set, operation 834d4bb536fSBarry Smith 835ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 836d4bb536fSBarry Smith @*/ 8377087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 838d4bb536fSBarry Smith { 839dfbe8321SBarry Smith PetscErrorCode ierr; 840ace3abfcSBarry Smith PetscBool flg; 841273d9f13SBarry Smith 8423a40ed3dSBarry Smith PetscFunctionBegin; 8430700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 844d4bb536fSBarry Smith if (op == MATOP_DESTROY) { 845251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 846273d9f13SBarry Smith if (flg) { 847d4bb536fSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 848c134de8dSSatish Balay *f = (void (*)(void))shell->destroy; 849c7fcc2eaSBarry Smith } else { 850c134de8dSSatish Balay *f = (void (*)(void))mat->ops->destroy; 851d4bb536fSBarry Smith } 852c7fcc2eaSBarry Smith } else if (op == MATOP_VIEW) { 853c134de8dSSatish Balay *f = (void (*)(void))mat->ops->view; 854c7fcc2eaSBarry Smith } else { 855c134de8dSSatish Balay *f = (((void (**)(void))mat->ops)[op]); 856d4bb536fSBarry Smith } 8573a40ed3dSBarry Smith PetscFunctionReturn(0); 858d4bb536fSBarry Smith } 859d4bb536fSBarry Smith 860