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 8c6db04a5SJed Brown #include <private/matimpl.h> /*I "petscmat.h" I*/ 9c6db04a5SJed Brown #include <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); 16ef66eb69SBarry Smith PetscScalar vscale,vshift; 178fe8eb27SJed Brown Vec dshift; 188fe8eb27SJed Brown Vec left,right; 198fe8eb27SJed Brown Vec dshift_owned,left_owned,right_owned; 208fe8eb27SJed Brown Vec left_work,right_work; 218fe8eb27SJed Brown PetscBool usingscaled; 2220563c6bSBarry Smith void *ctx; 2388cf3e7dSBarry Smith } Mat_Shell; 248fe8eb27SJed Brown /* 258fe8eb27SJed Brown The most general expression for the matrix is 268fe8eb27SJed Brown 278fe8eb27SJed Brown S = L (a A + B) R 288fe8eb27SJed Brown 298fe8eb27SJed Brown where 308fe8eb27SJed Brown A is the matrix defined by the user's function 318fe8eb27SJed Brown a is a scalar multiple 328fe8eb27SJed Brown L is left scaling 338fe8eb27SJed Brown R is right scaling 348fe8eb27SJed Brown B is a diagonal shift defined by 358fe8eb27SJed Brown diag(dshift) if the vector dshift is non-NULL 368fe8eb27SJed Brown vscale*identity otherwise 378fe8eb27SJed Brown 388fe8eb27SJed Brown The following identities apply: 398fe8eb27SJed Brown 408fe8eb27SJed Brown Scale by c: 418fe8eb27SJed Brown c [L (a A + B) R] = L [(a c) A + c B] R 428fe8eb27SJed Brown 438fe8eb27SJed Brown Shift by c: 448fe8eb27SJed Brown [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R 458fe8eb27SJed Brown 468fe8eb27SJed Brown Diagonal scaling is achieved by simply multiplying with existing L and R vectors 478fe8eb27SJed Brown 488fe8eb27SJed Brown In the data structure: 498fe8eb27SJed Brown 508fe8eb27SJed Brown vscale=1.0 means no special scaling will be applied 518fe8eb27SJed Brown dshift=NULL means a constant diagonal shift (fall back to vshift) 528fe8eb27SJed Brown vshift=0.0 means no constant diagonal shift, note that vshift is only used if dshift is NULL 538fe8eb27SJed Brown */ 548fe8eb27SJed Brown 558fe8eb27SJed Brown static PetscErrorCode MatMult_Shell(Mat,Vec,Vec); 568fe8eb27SJed Brown static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec); 578fe8eb27SJed Brown static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec); 588fe8eb27SJed Brown 598fe8eb27SJed Brown #undef __FUNCT__ 608fe8eb27SJed Brown #define __FUNCT__ "MatShellUseScaledMethods" 618fe8eb27SJed Brown static PetscErrorCode MatShellUseScaledMethods(Mat Y) 628fe8eb27SJed Brown { 638fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 648fe8eb27SJed Brown 658fe8eb27SJed Brown PetscFunctionBegin; 668fe8eb27SJed Brown if (shell->usingscaled) PetscFunctionReturn(0); 678fe8eb27SJed Brown shell->mult = Y->ops->mult; 688fe8eb27SJed Brown Y->ops->mult = MatMult_Shell; 698fe8eb27SJed Brown if (Y->ops->multtranspose) { 708fe8eb27SJed Brown shell->multtranspose = Y->ops->multtranspose; 718fe8eb27SJed Brown Y->ops->multtranspose = MatMultTranspose_Shell; 728fe8eb27SJed Brown } 738fe8eb27SJed Brown if (Y->ops->getdiagonal) { 748fe8eb27SJed Brown shell->getdiagonal = Y->ops->getdiagonal; 758fe8eb27SJed Brown Y->ops->getdiagonal = MatGetDiagonal_Shell; 768fe8eb27SJed Brown } 778fe8eb27SJed Brown shell->usingscaled = PETSC_TRUE; 788fe8eb27SJed Brown PetscFunctionReturn(0); 798fe8eb27SJed Brown } 808fe8eb27SJed Brown 818fe8eb27SJed Brown #undef __FUNCT__ 828fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleLeft" 838fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 848fe8eb27SJed Brown { 858fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 868fe8eb27SJed Brown PetscErrorCode ierr; 878fe8eb27SJed Brown 888fe8eb27SJed Brown PetscFunctionBegin; 898fe8eb27SJed Brown if (!shell->left) { 908fe8eb27SJed Brown *xx = x; 918fe8eb27SJed Brown } else { 928fe8eb27SJed Brown if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);} 938fe8eb27SJed Brown ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr); 948fe8eb27SJed Brown *xx = shell->left_work; 958fe8eb27SJed Brown } 968fe8eb27SJed Brown PetscFunctionReturn(0); 978fe8eb27SJed Brown } 988fe8eb27SJed Brown 998fe8eb27SJed Brown #undef __FUNCT__ 1008fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleRight" 1018fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 1028fe8eb27SJed Brown { 1038fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1048fe8eb27SJed Brown PetscErrorCode ierr; 1058fe8eb27SJed Brown 1068fe8eb27SJed Brown PetscFunctionBegin; 1078fe8eb27SJed Brown if (!shell->right) { 1088fe8eb27SJed Brown *xx = x; 1098fe8eb27SJed Brown } else { 1108fe8eb27SJed Brown if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);} 1118fe8eb27SJed Brown ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr); 1128fe8eb27SJed Brown *xx = shell->right_work; 1138fe8eb27SJed Brown } 1148fe8eb27SJed Brown PetscFunctionReturn(0); 1158fe8eb27SJed Brown } 1168fe8eb27SJed Brown 1178fe8eb27SJed Brown #undef __FUNCT__ 1188fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleLeft" 1198fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 1208fe8eb27SJed Brown { 1218fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1228fe8eb27SJed Brown PetscErrorCode ierr; 1238fe8eb27SJed Brown 1248fe8eb27SJed Brown PetscFunctionBegin; 1258fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);} 1268fe8eb27SJed Brown PetscFunctionReturn(0); 1278fe8eb27SJed Brown } 1288fe8eb27SJed Brown 1298fe8eb27SJed Brown #undef __FUNCT__ 1308fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleRight" 1318fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 1328fe8eb27SJed Brown { 1338fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1348fe8eb27SJed Brown PetscErrorCode ierr; 1358fe8eb27SJed Brown 1368fe8eb27SJed Brown PetscFunctionBegin; 1378fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);} 1388fe8eb27SJed Brown PetscFunctionReturn(0); 1398fe8eb27SJed Brown } 1408fe8eb27SJed Brown 1418fe8eb27SJed Brown #undef __FUNCT__ 1428fe8eb27SJed Brown #define __FUNCT__ "MatShellShiftAndScale" 1438fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 1448fe8eb27SJed Brown { 1458fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1468fe8eb27SJed Brown PetscErrorCode ierr; 1478fe8eb27SJed Brown 1488fe8eb27SJed Brown PetscFunctionBegin; 1498fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 1508fe8eb27SJed Brown PetscInt i,m; 1518fe8eb27SJed Brown const PetscScalar *x,*d; 1528fe8eb27SJed Brown PetscScalar *y; 1538fe8eb27SJed Brown ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr); 1548fe8eb27SJed Brown ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1558fe8eb27SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 1568fe8eb27SJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1578fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 1588fe8eb27SJed Brown ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1598fe8eb27SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 1608fe8eb27SJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 161*027c5fb5SJed Brown } else if (shell->vshift != 0) { 1628fe8eb27SJed Brown ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr); 163*027c5fb5SJed Brown } else if (shell->vscale != 1.0) { 164*027c5fb5SJed Brown ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr); 1658fe8eb27SJed Brown } 1668fe8eb27SJed Brown PetscFunctionReturn(0); 1678fe8eb27SJed Brown } 168e51e0e81SBarry Smith 169711e205bSSatish Balay #undef __FUNCT__ 170711e205bSSatish Balay #define __FUNCT__ "MatShellGetContext" 17135d520bfSBarry Smith /*@C 172a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 173b4fd4287SBarry Smith 174c7fcc2eaSBarry Smith Not Collective 175c7fcc2eaSBarry Smith 176b4fd4287SBarry Smith Input Parameter: 177b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 178b4fd4287SBarry Smith 179b4fd4287SBarry Smith Output Parameter: 180b4fd4287SBarry Smith . ctx - the user provided context 181b4fd4287SBarry Smith 18215091d37SBarry Smith Level: advanced 18315091d37SBarry Smith 184a62d957aSLois Curfman McInnes Notes: 185a62d957aSLois Curfman McInnes This routine is intended for use within various shell matrix routines, 186a62d957aSLois Curfman McInnes as set with MatShellSetOperation(). 187a62d957aSLois Curfman McInnes 188a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context 189a62d957aSLois Curfman McInnes 190ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 1910bc0a719SBarry Smith @*/ 1928fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 193b4fd4287SBarry Smith { 194dfbe8321SBarry Smith PetscErrorCode ierr; 195ace3abfcSBarry Smith PetscBool flg; 196273d9f13SBarry Smith 1973a40ed3dSBarry Smith PetscFunctionBegin; 1980700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1994482741eSBarry Smith PetscValidPointer(ctx,2); 200273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 2018fe8eb27SJed Brown if (!flg) *(void**)ctx = 0; 2028fe8eb27SJed Brown else *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 2033a40ed3dSBarry Smith PetscFunctionReturn(0); 204b4fd4287SBarry Smith } 205b4fd4287SBarry Smith 206711e205bSSatish Balay #undef __FUNCT__ 207711e205bSSatish Balay #define __FUNCT__ "MatDestroy_Shell" 208dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 209e51e0e81SBarry Smith { 210dfbe8321SBarry Smith PetscErrorCode ierr; 211bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 212ed3cc1f0SBarry Smith 2133a40ed3dSBarry Smith PetscFunctionBegin; 214bf0cc555SLisandro Dalcin if (shell->destroy) { 215bf0cc555SLisandro Dalcin ierr = (*shell->destroy)(mat);CHKERRQ(ierr); 216bf0cc555SLisandro Dalcin } 2178fe8eb27SJed Brown ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr); 2188fe8eb27SJed Brown ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr); 2198fe8eb27SJed Brown ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr); 2208fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 2218fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 222bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 2233a40ed3dSBarry Smith PetscFunctionReturn(0); 224e51e0e81SBarry Smith } 225e51e0e81SBarry Smith 226711e205bSSatish Balay #undef __FUNCT__ 227ef66eb69SBarry Smith #define __FUNCT__ "MatMult_Shell" 228dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 229ef66eb69SBarry Smith { 230ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 231dfbe8321SBarry Smith PetscErrorCode ierr; 2328fe8eb27SJed Brown Vec xx; 233ef66eb69SBarry Smith 234ef66eb69SBarry Smith PetscFunctionBegin; 2358fe8eb27SJed Brown ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr); 2368fe8eb27SJed Brown ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr); 2378fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2388fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 239ef66eb69SBarry Smith PetscFunctionReturn(0); 240ef66eb69SBarry Smith } 241ef66eb69SBarry Smith 242ef66eb69SBarry Smith #undef __FUNCT__ 24318be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell" 24418be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 24518be62a5SSatish Balay { 24618be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 24718be62a5SSatish Balay PetscErrorCode ierr; 2488fe8eb27SJed Brown Vec xx; 24918be62a5SSatish Balay 25018be62a5SSatish Balay PetscFunctionBegin; 2518fe8eb27SJed Brown ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr); 2528fe8eb27SJed Brown ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr); 2538fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2548fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 25518be62a5SSatish Balay PetscFunctionReturn(0); 25618be62a5SSatish Balay } 25718be62a5SSatish Balay 25818be62a5SSatish Balay #undef __FUNCT__ 25918be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell" 26018be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 26118be62a5SSatish Balay { 26218be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 26318be62a5SSatish Balay PetscErrorCode ierr; 26418be62a5SSatish Balay 26518be62a5SSatish Balay PetscFunctionBegin; 26618be62a5SSatish Balay ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr); 26718be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 2688fe8eb27SJed Brown if (shell->dshift) { 2698fe8eb27SJed Brown ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr); 2708fe8eb27SJed Brown } else { 27118be62a5SSatish Balay ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 27218be62a5SSatish Balay } 2738fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 2748fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 27518be62a5SSatish Balay PetscFunctionReturn(0); 27618be62a5SSatish Balay } 27718be62a5SSatish Balay 27818be62a5SSatish Balay #undef __FUNCT__ 279ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell" 280f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 281ef66eb69SBarry Smith { 282ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 2838fe8eb27SJed Brown PetscErrorCode ierr; 284b24ad042SBarry Smith 285ef66eb69SBarry Smith PetscFunctionBegin; 2868fe8eb27SJed Brown if (shell->left || shell->right || shell->dshift) { 2878fe8eb27SJed Brown if (!shell->dshift) { 2888fe8eb27SJed Brown if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);} 2898fe8eb27SJed Brown shell->dshift = shell->dshift_owned; 2908fe8eb27SJed Brown ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr); 2918fe8eb27SJed Brown } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);} 2928fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 2938fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 2948fe8eb27SJed Brown } else shell->vshift += a; 2958fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 296ef66eb69SBarry Smith PetscFunctionReturn(0); 297ef66eb69SBarry Smith } 298ef66eb69SBarry Smith 299ef66eb69SBarry Smith #undef __FUNCT__ 300ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell" 301f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 302ef66eb69SBarry Smith { 303ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3048fe8eb27SJed Brown PetscErrorCode ierr; 305b24ad042SBarry Smith 306ef66eb69SBarry Smith PetscFunctionBegin; 307f4df32b1SMatthew Knepley shell->vscale *= a; 3088fe8eb27SJed Brown if (shell->dshift) {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);} 3098fe8eb27SJed Brown else shell->vshift *= a; 3108fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 3118fe8eb27SJed Brown PetscFunctionReturn(0); 31218be62a5SSatish Balay } 3138fe8eb27SJed Brown 3148fe8eb27SJed Brown #undef __FUNCT__ 3158fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell" 3168fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 3178fe8eb27SJed Brown { 3188fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 3198fe8eb27SJed Brown PetscErrorCode ierr; 3208fe8eb27SJed Brown 3218fe8eb27SJed Brown PetscFunctionBegin; 3228fe8eb27SJed Brown if (left) { 3238fe8eb27SJed Brown if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);} 3248fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);} 3258fe8eb27SJed Brown else { 3268fe8eb27SJed Brown shell->left = shell->left_owned; 3278fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 32818be62a5SSatish Balay } 329ef66eb69SBarry Smith } 3308fe8eb27SJed Brown if (right) { 3318fe8eb27SJed Brown if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);} 3328fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);} 3338fe8eb27SJed Brown else { 3348fe8eb27SJed Brown shell->right = shell->right_owned; 3358fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 3368fe8eb27SJed Brown } 3378fe8eb27SJed Brown } 3388fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 339ef66eb69SBarry Smith PetscFunctionReturn(0); 340ef66eb69SBarry Smith } 341ef66eb69SBarry Smith 342ef66eb69SBarry Smith #undef __FUNCT__ 343ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell" 344dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 345ef66eb69SBarry Smith { 346ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 347ef66eb69SBarry Smith 348ef66eb69SBarry Smith PetscFunctionBegin; 3498fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 350ef66eb69SBarry Smith shell->vshift = 0.0; 351ef66eb69SBarry Smith shell->vscale = 1.0; 3528fe8eb27SJed Brown shell->dshift = PETSC_NULL; 3538fe8eb27SJed Brown shell->left = PETSC_NULL; 3548fe8eb27SJed Brown shell->right = PETSC_NULL; 3557fb7d96aSJed Brown if (shell->mult) { 356ef66eb69SBarry Smith Y->ops->mult = shell->mult; 3577fb7d96aSJed Brown shell->mult = PETSC_NULL; 3587fb7d96aSJed Brown } 3597fb7d96aSJed Brown if (shell->multtranspose) { 36018be62a5SSatish Balay Y->ops->multtranspose = shell->multtranspose; 3617fb7d96aSJed Brown shell->multtranspose = PETSC_NULL; 3627fb7d96aSJed Brown } 3637fb7d96aSJed Brown if (shell->getdiagonal) { 36418be62a5SSatish Balay Y->ops->getdiagonal = shell->getdiagonal; 3657fb7d96aSJed Brown shell->getdiagonal = PETSC_NULL; 3667fb7d96aSJed Brown } 3677fb7d96aSJed Brown shell->usingscaled = PETSC_FALSE; 368ef66eb69SBarry Smith } 369ef66eb69SBarry Smith PetscFunctionReturn(0); 370ef66eb69SBarry Smith } 371ef66eb69SBarry Smith 37209573ac7SBarry Smith extern PetscErrorCode MatConvert_Shell(Mat, const MatType,MatReuse,Mat*); 373b951964fSBarry Smith 374521d7252SBarry Smith #undef __FUNCT__ 375521d7252SBarry Smith #define __FUNCT__ "MatSetBlockSize_Shell" 376521d7252SBarry Smith PetscErrorCode MatSetBlockSize_Shell(Mat A,PetscInt bs) 377521d7252SBarry Smith { 37841c166b1SJed Brown PetscErrorCode ierr; 37941c166b1SJed Brown 380521d7252SBarry Smith PetscFunctionBegin; 38141c166b1SJed Brown ierr = PetscLayoutSetBlockSize(A->rmap,bs);CHKERRQ(ierr); 38241c166b1SJed Brown ierr = PetscLayoutSetBlockSize(A->cmap,bs);CHKERRQ(ierr); 383521d7252SBarry Smith PetscFunctionReturn(0); 384521d7252SBarry Smith } 385521d7252SBarry Smith 38609dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 38720563c6bSBarry Smith 0, 38820563c6bSBarry Smith 0, 38920563c6bSBarry Smith 0, 39097304618SKris Buschelman /* 4*/ 0, 39120563c6bSBarry Smith 0, 392b951964fSBarry Smith 0, 393b951964fSBarry Smith 0, 394b951964fSBarry Smith 0, 395b951964fSBarry Smith 0, 39697304618SKris Buschelman /*10*/ 0, 397b951964fSBarry Smith 0, 398b951964fSBarry Smith 0, 399b951964fSBarry Smith 0, 400b951964fSBarry Smith 0, 40197304618SKris Buschelman /*15*/ 0, 402b951964fSBarry Smith 0, 403b951964fSBarry Smith 0, 4048fe8eb27SJed Brown MatDiagonalScale_Shell, 405b951964fSBarry Smith 0, 40697304618SKris Buschelman /*20*/ 0, 407ef66eb69SBarry Smith MatAssemblyEnd_Shell, 408b951964fSBarry Smith 0, 409b951964fSBarry Smith 0, 410d519adbfSMatthew Knepley /*24*/ 0, 411b951964fSBarry Smith 0, 412b951964fSBarry Smith 0, 413b951964fSBarry Smith 0, 414b951964fSBarry Smith 0, 415d519adbfSMatthew Knepley /*29*/ 0, 416b951964fSBarry Smith 0, 417273d9f13SBarry Smith 0, 418b951964fSBarry Smith 0, 419b951964fSBarry Smith 0, 420d519adbfSMatthew Knepley /*34*/ 0, 421b951964fSBarry Smith 0, 422b951964fSBarry Smith 0, 42309dc0095SBarry Smith 0, 42409dc0095SBarry Smith 0, 425d519adbfSMatthew Knepley /*39*/ 0, 42609dc0095SBarry Smith 0, 42709dc0095SBarry Smith 0, 42809dc0095SBarry Smith 0, 42909dc0095SBarry Smith 0, 430d519adbfSMatthew Knepley /*44*/ 0, 431ef66eb69SBarry Smith MatScale_Shell, 432ef66eb69SBarry Smith MatShift_Shell, 43309dc0095SBarry Smith 0, 43409dc0095SBarry Smith 0, 435d519adbfSMatthew Knepley /*49*/ MatSetBlockSize_Shell, 43609dc0095SBarry Smith 0, 43709dc0095SBarry Smith 0, 43809dc0095SBarry Smith 0, 43909dc0095SBarry Smith 0, 440d519adbfSMatthew Knepley /*54*/ 0, 44109dc0095SBarry Smith 0, 44209dc0095SBarry Smith 0, 44309dc0095SBarry Smith 0, 44409dc0095SBarry Smith 0, 445d519adbfSMatthew Knepley /*59*/ 0, 446b9b97703SBarry Smith MatDestroy_Shell, 44709dc0095SBarry Smith 0, 448357abbc8SBarry Smith 0, 449273d9f13SBarry Smith 0, 450d519adbfSMatthew Knepley /*64*/ 0, 451273d9f13SBarry Smith 0, 452273d9f13SBarry Smith 0, 453273d9f13SBarry Smith 0, 454273d9f13SBarry Smith 0, 455d519adbfSMatthew Knepley /*69*/ 0, 456273d9f13SBarry Smith 0, 457c87e5d42SMatthew Knepley MatConvert_Shell, 458273d9f13SBarry Smith 0, 45997304618SKris Buschelman 0, 460d519adbfSMatthew Knepley /*74*/ 0, 46197304618SKris Buschelman 0, 46297304618SKris Buschelman 0, 46397304618SKris Buschelman 0, 46497304618SKris Buschelman 0, 465d519adbfSMatthew Knepley /*79*/ 0, 46697304618SKris Buschelman 0, 46797304618SKris Buschelman 0, 46897304618SKris Buschelman 0, 469865e5f61SKris Buschelman 0, 470d519adbfSMatthew Knepley /*84*/ 0, 471865e5f61SKris Buschelman 0, 472865e5f61SKris Buschelman 0, 473865e5f61SKris Buschelman 0, 474865e5f61SKris Buschelman 0, 475d519adbfSMatthew Knepley /*89*/ 0, 476865e5f61SKris Buschelman 0, 477865e5f61SKris Buschelman 0, 478865e5f61SKris Buschelman 0, 479865e5f61SKris Buschelman 0, 480d519adbfSMatthew Knepley /*94*/ 0, 481865e5f61SKris Buschelman 0, 482865e5f61SKris Buschelman 0, 483865e5f61SKris Buschelman 0}; 484273d9f13SBarry Smith 4850bad9183SKris Buschelman /*MC 486fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 4870bad9183SKris Buschelman 4880bad9183SKris Buschelman Level: advanced 4890bad9183SKris Buschelman 4900bad9183SKris Buschelman .seealso: MatCreateShell 4910bad9183SKris Buschelman M*/ 4920bad9183SKris Buschelman 493273d9f13SBarry Smith EXTERN_C_BEGIN 494711e205bSSatish Balay #undef __FUNCT__ 495711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell" 4967087cfbeSBarry Smith PetscErrorCode MatCreate_Shell(Mat A) 497273d9f13SBarry Smith { 498273d9f13SBarry Smith Mat_Shell *b; 499dfbe8321SBarry Smith PetscErrorCode ierr; 500273d9f13SBarry Smith 501273d9f13SBarry Smith PetscFunctionBegin; 502273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 503273d9f13SBarry Smith 50438f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_Shell,&b);CHKERRQ(ierr); 505273d9f13SBarry Smith A->data = (void*)b; 506273d9f13SBarry Smith 50726283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 50826283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 50926283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 51026283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 511273d9f13SBarry Smith 512273d9f13SBarry Smith b->ctx = 0; 513ef66eb69SBarry Smith b->vshift = 0.0; 514ef66eb69SBarry Smith b->vscale = 1.0; 515ef66eb69SBarry Smith b->mult = 0; 51618be62a5SSatish Balay b->multtranspose = 0; 51718be62a5SSatish Balay b->getdiagonal = 0; 518273d9f13SBarry Smith A->assembled = PETSC_TRUE; 519210f0121SBarry Smith A->preallocated = PETSC_FALSE; 52017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 521273d9f13SBarry Smith PetscFunctionReturn(0); 522273d9f13SBarry Smith } 523273d9f13SBarry Smith EXTERN_C_END 524e51e0e81SBarry Smith 525711e205bSSatish Balay #undef __FUNCT__ 526711e205bSSatish Balay #define __FUNCT__ "MatCreateShell" 5274b828684SBarry Smith /*@C 528052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 529ff756334SLois Curfman McInnes private data storage format. 530e51e0e81SBarry Smith 531c7fcc2eaSBarry Smith Collective on MPI_Comm 532c7fcc2eaSBarry Smith 533e51e0e81SBarry Smith Input Parameters: 534c7fcc2eaSBarry Smith + comm - MPI communicator 535c7fcc2eaSBarry Smith . m - number of local rows (must be given) 536c7fcc2eaSBarry Smith . n - number of local columns (must be given) 537c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 538c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 539c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 540e51e0e81SBarry Smith 541ff756334SLois Curfman McInnes Output Parameter: 54244cd7ae7SLois Curfman McInnes . A - the matrix 543e51e0e81SBarry Smith 544ff2fd236SBarry Smith Level: advanced 545ff2fd236SBarry Smith 546f39d1f56SLois Curfman McInnes Usage: 5477b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 548f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 549c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 550f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 551f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 552f39d1f56SLois Curfman McInnes 553ff756334SLois Curfman McInnes Notes: 554ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 555ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 556ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 557e51e0e81SBarry Smith 558f38a8d56SBarry Smith Fortran Notes: The context can only be an integer or a PetscObject 559f38a8d56SBarry Smith unfortunately it cannot be a Fortran array or derived type. 560f38a8d56SBarry Smith 561f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 562f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 563645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 564645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 565645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 566645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 567645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 568645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 569645985a0SLois Curfman McInnes For example, 570f39d1f56SLois Curfman McInnes 571f39d1f56SLois Curfman McInnes $ 572f39d1f56SLois Curfman McInnes $ Vec x, y 5737b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 574645985a0SLois Curfman McInnes $ Mat A 575f39d1f56SLois Curfman McInnes $ 576c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 577c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 578f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 579c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 580c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 581c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 582645985a0SLois Curfman McInnes $ MatMult(A,x,y); 583645985a0SLois Curfman McInnes $ MatDestroy(A); 584f39d1f56SLois Curfman McInnes $ VecDestroy(y); VecDestroy(x); 585645985a0SLois Curfman McInnes $ 586e51e0e81SBarry Smith 5870b627109SLois Curfman McInnes .keywords: matrix, shell, create 5880b627109SLois Curfman McInnes 589ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext() 590e51e0e81SBarry Smith @*/ 5917087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 592e51e0e81SBarry Smith { 593dfbe8321SBarry Smith PetscErrorCode ierr; 594ed3cc1f0SBarry Smith 5953a40ed3dSBarry Smith PetscFunctionBegin; 596f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 597f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 598899cda47SBarry Smith 599273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 600273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 601273d9f13SBarry Smith PetscFunctionReturn(0); 602c7fcc2eaSBarry Smith } 603c7fcc2eaSBarry Smith 604711e205bSSatish Balay #undef __FUNCT__ 605711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext" 606c6866cfdSSatish Balay /*@ 607273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 608c7fcc2eaSBarry Smith 6093f9fe445SBarry Smith Logically Collective on Mat 610c7fcc2eaSBarry Smith 611273d9f13SBarry Smith Input Parameters: 612273d9f13SBarry Smith + mat - the shell matrix 613273d9f13SBarry Smith - ctx - the context 614273d9f13SBarry Smith 615273d9f13SBarry Smith Level: advanced 616273d9f13SBarry Smith 617f38a8d56SBarry Smith Fortran Notes: The context can only be an integer or a PetscObject 618f38a8d56SBarry Smith unfortunately it cannot be a Fortran array or derived type. 619273d9f13SBarry Smith 620273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 6210bc0a719SBarry Smith @*/ 6227087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 623273d9f13SBarry Smith { 624273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 625dfbe8321SBarry Smith PetscErrorCode ierr; 626ace3abfcSBarry Smith PetscBool flg; 627273d9f13SBarry Smith 628273d9f13SBarry Smith PetscFunctionBegin; 6290700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 630273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 631273d9f13SBarry Smith if (flg) { 632273d9f13SBarry Smith shell->ctx = ctx; 633273d9f13SBarry Smith } 6343a40ed3dSBarry Smith PetscFunctionReturn(0); 635e51e0e81SBarry Smith } 636e51e0e81SBarry Smith 637711e205bSSatish Balay #undef __FUNCT__ 638711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation" 639c16cb8f2SBarry Smith /*@C 6403a3eedf2SBarry Smith MatShellSetOperation - Allows user to set a matrix operation for 6413a3eedf2SBarry Smith a shell matrix. 642e51e0e81SBarry Smith 6433f9fe445SBarry Smith Logically Collective on Mat 644fee21e36SBarry Smith 645c7fcc2eaSBarry Smith Input Parameters: 646c7fcc2eaSBarry Smith + mat - the shell matrix 647c7fcc2eaSBarry Smith . op - the name of the operation 648c7fcc2eaSBarry Smith - f - the function that provides the operation. 649c7fcc2eaSBarry Smith 65015091d37SBarry Smith Level: advanced 65115091d37SBarry Smith 652fae171e0SBarry Smith Usage: 653e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 654f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 655c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 6560b627109SLois Curfman McInnes 657a62d957aSLois Curfman McInnes Notes: 658e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 6591c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 660a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 6611c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 662a62d957aSLois Curfman McInnes 663a62d957aSLois Curfman McInnes All user-provided functions should have the same calling 664deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 665deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 666deebb3c3SLois Curfman McInnes routines, e.g., 667a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 668a62d957aSLois Curfman McInnes 6694aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 6704aa34b0aSBarry Smith nonzero on failure. 6714aa34b0aSBarry Smith 672a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 673a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 674a62d957aSLois Curfman McInnes set by MatCreateShell(). 675a62d957aSLois Curfman McInnes 676501d9185SBarry Smith Fortran Notes: For MatGetVecs() the user code should check if the input left or right matrix is -1 and in that case not 677501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 678501d9185SBarry Smith 679a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation 680a62d957aSLois Curfman McInnes 681ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext() 682e51e0e81SBarry Smith @*/ 6837087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 684e51e0e81SBarry Smith { 685dfbe8321SBarry Smith PetscErrorCode ierr; 686ace3abfcSBarry Smith PetscBool flg; 687273d9f13SBarry Smith 6883a40ed3dSBarry Smith PetscFunctionBegin; 6890700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6901c1c02c0SLois Curfman McInnes if (op == MATOP_DESTROY) { 691273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 692273d9f13SBarry Smith if (flg) { 693a62d957aSLois Curfman McInnes Mat_Shell *shell = (Mat_Shell*)mat->data; 6946849ba73SBarry Smith shell->destroy = (PetscErrorCode (*)(Mat)) f; 6956849ba73SBarry Smith } else mat->ops->destroy = (PetscErrorCode (*)(Mat)) f; 696a62d957aSLois Curfman McInnes } 6976849ba73SBarry Smith else if (op == MATOP_VIEW) mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer)) f; 698c134de8dSSatish Balay else (((void(**)(void))mat->ops)[op]) = f; 699a62d957aSLois Curfman McInnes 7003a40ed3dSBarry Smith PetscFunctionReturn(0); 701e51e0e81SBarry Smith } 702f0479e8cSBarry Smith 703711e205bSSatish Balay #undef __FUNCT__ 704711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation" 705d4bb536fSBarry Smith /*@C 706d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 707d4bb536fSBarry Smith 708c7fcc2eaSBarry Smith Not Collective 709c7fcc2eaSBarry Smith 710d4bb536fSBarry Smith Input Parameters: 711c7fcc2eaSBarry Smith + mat - the shell matrix 712c7fcc2eaSBarry Smith - op - the name of the operation 713d4bb536fSBarry Smith 714d4bb536fSBarry Smith Output Parameter: 715d4bb536fSBarry Smith . f - the function that provides the operation. 716d4bb536fSBarry Smith 71715091d37SBarry Smith Level: advanced 71815091d37SBarry Smith 719d4bb536fSBarry Smith Notes: 720e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 721d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 722d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 723d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 724d4bb536fSBarry Smith 725d4bb536fSBarry Smith All user-provided functions have the same calling 726d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 727d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 728d4bb536fSBarry Smith routines, e.g., 729d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 730d4bb536fSBarry Smith 731d4bb536fSBarry Smith Within each user-defined routine, the user should call 732d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 733d4bb536fSBarry Smith set by MatCreateShell(). 734d4bb536fSBarry Smith 735d4bb536fSBarry Smith .keywords: matrix, shell, set, operation 736d4bb536fSBarry Smith 737ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 738d4bb536fSBarry Smith @*/ 7397087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 740d4bb536fSBarry Smith { 741dfbe8321SBarry Smith PetscErrorCode ierr; 742ace3abfcSBarry Smith PetscBool flg; 743273d9f13SBarry Smith 7443a40ed3dSBarry Smith PetscFunctionBegin; 7450700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 746d4bb536fSBarry Smith if (op == MATOP_DESTROY) { 747273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 748273d9f13SBarry Smith if (flg) { 749d4bb536fSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 750c134de8dSSatish Balay *f = (void(*)(void))shell->destroy; 751c7fcc2eaSBarry Smith } else { 752c134de8dSSatish Balay *f = (void(*)(void))mat->ops->destroy; 753d4bb536fSBarry Smith } 754c7fcc2eaSBarry Smith } else if (op == MATOP_VIEW) { 755c134de8dSSatish Balay *f = (void(*)(void))mat->ops->view; 756c7fcc2eaSBarry Smith } else { 757c134de8dSSatish Balay *f = (((void(**)(void))mat->ops)[op]); 758d4bb536fSBarry Smith } 759d4bb536fSBarry Smith 7603a40ed3dSBarry Smith PetscFunctionReturn(0); 761d4bb536fSBarry Smith } 762d4bb536fSBarry Smith 763