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; 17*8fe8eb27SJed Brown Vec dshift; 18*8fe8eb27SJed Brown Vec left,right; 19*8fe8eb27SJed Brown Vec dshift_owned,left_owned,right_owned; 20*8fe8eb27SJed Brown Vec left_work,right_work; 21*8fe8eb27SJed Brown PetscBool usingscaled; 2220563c6bSBarry Smith void *ctx; 2388cf3e7dSBarry Smith } Mat_Shell; 24*8fe8eb27SJed Brown /* 25*8fe8eb27SJed Brown The most general expression for the matrix is 26*8fe8eb27SJed Brown 27*8fe8eb27SJed Brown S = L (a A + B) R 28*8fe8eb27SJed Brown 29*8fe8eb27SJed Brown where 30*8fe8eb27SJed Brown A is the matrix defined by the user's function 31*8fe8eb27SJed Brown a is a scalar multiple 32*8fe8eb27SJed Brown L is left scaling 33*8fe8eb27SJed Brown R is right scaling 34*8fe8eb27SJed Brown B is a diagonal shift defined by 35*8fe8eb27SJed Brown diag(dshift) if the vector dshift is non-NULL 36*8fe8eb27SJed Brown vscale*identity otherwise 37*8fe8eb27SJed Brown 38*8fe8eb27SJed Brown The following identities apply: 39*8fe8eb27SJed Brown 40*8fe8eb27SJed Brown Scale by c: 41*8fe8eb27SJed Brown c [L (a A + B) R] = L [(a c) A + c B] R 42*8fe8eb27SJed Brown 43*8fe8eb27SJed Brown Shift by c: 44*8fe8eb27SJed Brown [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R 45*8fe8eb27SJed Brown 46*8fe8eb27SJed Brown Diagonal scaling is achieved by simply multiplying with existing L and R vectors 47*8fe8eb27SJed Brown 48*8fe8eb27SJed Brown In the data structure: 49*8fe8eb27SJed Brown 50*8fe8eb27SJed Brown vscale=1.0 means no special scaling will be applied 51*8fe8eb27SJed Brown dshift=NULL means a constant diagonal shift (fall back to vshift) 52*8fe8eb27SJed Brown vshift=0.0 means no constant diagonal shift, note that vshift is only used if dshift is NULL 53*8fe8eb27SJed Brown */ 54*8fe8eb27SJed Brown 55*8fe8eb27SJed Brown static PetscErrorCode MatMult_Shell(Mat,Vec,Vec); 56*8fe8eb27SJed Brown static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec); 57*8fe8eb27SJed Brown static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec); 58*8fe8eb27SJed Brown 59*8fe8eb27SJed Brown #undef __FUNCT__ 60*8fe8eb27SJed Brown #define __FUNCT__ "MatShellUseScaledMethods" 61*8fe8eb27SJed Brown static PetscErrorCode MatShellUseScaledMethods(Mat Y) 62*8fe8eb27SJed Brown { 63*8fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 64*8fe8eb27SJed Brown 65*8fe8eb27SJed Brown PetscFunctionBegin; 66*8fe8eb27SJed Brown if (shell->usingscaled) PetscFunctionReturn(0); 67*8fe8eb27SJed Brown shell->mult = Y->ops->mult; 68*8fe8eb27SJed Brown Y->ops->mult = MatMult_Shell; 69*8fe8eb27SJed Brown if (Y->ops->multtranspose) { 70*8fe8eb27SJed Brown shell->multtranspose = Y->ops->multtranspose; 71*8fe8eb27SJed Brown Y->ops->multtranspose = MatMultTranspose_Shell; 72*8fe8eb27SJed Brown } 73*8fe8eb27SJed Brown if (Y->ops->getdiagonal) { 74*8fe8eb27SJed Brown shell->getdiagonal = Y->ops->getdiagonal; 75*8fe8eb27SJed Brown Y->ops->getdiagonal = MatGetDiagonal_Shell; 76*8fe8eb27SJed Brown } 77*8fe8eb27SJed Brown shell->usingscaled = PETSC_TRUE; 78*8fe8eb27SJed Brown PetscFunctionReturn(0); 79*8fe8eb27SJed Brown } 80*8fe8eb27SJed Brown 81*8fe8eb27SJed Brown #undef __FUNCT__ 82*8fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleLeft" 83*8fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 84*8fe8eb27SJed Brown { 85*8fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 86*8fe8eb27SJed Brown PetscErrorCode ierr; 87*8fe8eb27SJed Brown 88*8fe8eb27SJed Brown PetscFunctionBegin; 89*8fe8eb27SJed Brown if (!shell->left) { 90*8fe8eb27SJed Brown *xx = x; 91*8fe8eb27SJed Brown } else { 92*8fe8eb27SJed Brown if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);} 93*8fe8eb27SJed Brown ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr); 94*8fe8eb27SJed Brown *xx = shell->left_work; 95*8fe8eb27SJed Brown } 96*8fe8eb27SJed Brown PetscFunctionReturn(0); 97*8fe8eb27SJed Brown } 98*8fe8eb27SJed Brown 99*8fe8eb27SJed Brown #undef __FUNCT__ 100*8fe8eb27SJed Brown #define __FUNCT__ "MatShellPreScaleRight" 101*8fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 102*8fe8eb27SJed Brown { 103*8fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 104*8fe8eb27SJed Brown PetscErrorCode ierr; 105*8fe8eb27SJed Brown 106*8fe8eb27SJed Brown PetscFunctionBegin; 107*8fe8eb27SJed Brown if (!shell->right) { 108*8fe8eb27SJed Brown *xx = x; 109*8fe8eb27SJed Brown } else { 110*8fe8eb27SJed Brown if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);} 111*8fe8eb27SJed Brown ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr); 112*8fe8eb27SJed Brown *xx = shell->right_work; 113*8fe8eb27SJed Brown } 114*8fe8eb27SJed Brown PetscFunctionReturn(0); 115*8fe8eb27SJed Brown } 116*8fe8eb27SJed Brown 117*8fe8eb27SJed Brown #undef __FUNCT__ 118*8fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleLeft" 119*8fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 120*8fe8eb27SJed Brown { 121*8fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 122*8fe8eb27SJed Brown PetscErrorCode ierr; 123*8fe8eb27SJed Brown 124*8fe8eb27SJed Brown PetscFunctionBegin; 125*8fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);} 126*8fe8eb27SJed Brown PetscFunctionReturn(0); 127*8fe8eb27SJed Brown } 128*8fe8eb27SJed Brown 129*8fe8eb27SJed Brown #undef __FUNCT__ 130*8fe8eb27SJed Brown #define __FUNCT__ "MatShellPostScaleRight" 131*8fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 132*8fe8eb27SJed Brown { 133*8fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 134*8fe8eb27SJed Brown PetscErrorCode ierr; 135*8fe8eb27SJed Brown 136*8fe8eb27SJed Brown PetscFunctionBegin; 137*8fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);} 138*8fe8eb27SJed Brown PetscFunctionReturn(0); 139*8fe8eb27SJed Brown } 140*8fe8eb27SJed Brown 141*8fe8eb27SJed Brown #undef __FUNCT__ 142*8fe8eb27SJed Brown #define __FUNCT__ "MatShellShiftAndScale" 143*8fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 144*8fe8eb27SJed Brown { 145*8fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 146*8fe8eb27SJed Brown PetscErrorCode ierr; 147*8fe8eb27SJed Brown 148*8fe8eb27SJed Brown PetscFunctionBegin; 149*8fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 150*8fe8eb27SJed Brown PetscInt i,m; 151*8fe8eb27SJed Brown const PetscScalar *x,*d; 152*8fe8eb27SJed Brown PetscScalar *y; 153*8fe8eb27SJed Brown ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr); 154*8fe8eb27SJed Brown ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr); 155*8fe8eb27SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 156*8fe8eb27SJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 157*8fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 158*8fe8eb27SJed Brown ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr); 159*8fe8eb27SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 160*8fe8eb27SJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 161*8fe8eb27SJed Brown } else { 162*8fe8eb27SJed Brown ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr); 163*8fe8eb27SJed Brown } 164*8fe8eb27SJed Brown PetscFunctionReturn(0); 165*8fe8eb27SJed Brown } 166e51e0e81SBarry Smith 167711e205bSSatish Balay #undef __FUNCT__ 168711e205bSSatish Balay #define __FUNCT__ "MatShellGetContext" 16935d520bfSBarry Smith /*@C 170a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 171b4fd4287SBarry Smith 172c7fcc2eaSBarry Smith Not Collective 173c7fcc2eaSBarry Smith 174b4fd4287SBarry Smith Input Parameter: 175b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 176b4fd4287SBarry Smith 177b4fd4287SBarry Smith Output Parameter: 178b4fd4287SBarry Smith . ctx - the user provided context 179b4fd4287SBarry Smith 18015091d37SBarry Smith Level: advanced 18115091d37SBarry Smith 182a62d957aSLois Curfman McInnes Notes: 183a62d957aSLois Curfman McInnes This routine is intended for use within various shell matrix routines, 184a62d957aSLois Curfman McInnes as set with MatShellSetOperation(). 185a62d957aSLois Curfman McInnes 186a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context 187a62d957aSLois Curfman McInnes 188ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 1890bc0a719SBarry Smith @*/ 190*8fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 191b4fd4287SBarry Smith { 192dfbe8321SBarry Smith PetscErrorCode ierr; 193ace3abfcSBarry Smith PetscBool flg; 194273d9f13SBarry Smith 1953a40ed3dSBarry Smith PetscFunctionBegin; 1960700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1974482741eSBarry Smith PetscValidPointer(ctx,2); 198273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 199*8fe8eb27SJed Brown if (!flg) *(void**)ctx = 0; 200*8fe8eb27SJed Brown else *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 2013a40ed3dSBarry Smith PetscFunctionReturn(0); 202b4fd4287SBarry Smith } 203b4fd4287SBarry Smith 204711e205bSSatish Balay #undef __FUNCT__ 205711e205bSSatish Balay #define __FUNCT__ "MatDestroy_Shell" 206dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 207e51e0e81SBarry Smith { 208dfbe8321SBarry Smith PetscErrorCode ierr; 209bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 210ed3cc1f0SBarry Smith 2113a40ed3dSBarry Smith PetscFunctionBegin; 212bf0cc555SLisandro Dalcin if (shell->destroy) { 213bf0cc555SLisandro Dalcin ierr = (*shell->destroy)(mat);CHKERRQ(ierr); 214bf0cc555SLisandro Dalcin } 215*8fe8eb27SJed Brown ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr); 216*8fe8eb27SJed Brown ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr); 217*8fe8eb27SJed Brown ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr); 218*8fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 219*8fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 220bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 2213a40ed3dSBarry Smith PetscFunctionReturn(0); 222e51e0e81SBarry Smith } 223e51e0e81SBarry Smith 224711e205bSSatish Balay #undef __FUNCT__ 225ef66eb69SBarry Smith #define __FUNCT__ "MatMult_Shell" 226dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 227ef66eb69SBarry Smith { 228ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 229dfbe8321SBarry Smith PetscErrorCode ierr; 230*8fe8eb27SJed Brown Vec xx; 231ef66eb69SBarry Smith 232ef66eb69SBarry Smith PetscFunctionBegin; 233*8fe8eb27SJed Brown ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr); 234*8fe8eb27SJed Brown ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr); 235*8fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 236*8fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 237ef66eb69SBarry Smith PetscFunctionReturn(0); 238ef66eb69SBarry Smith } 239ef66eb69SBarry Smith 240ef66eb69SBarry Smith #undef __FUNCT__ 24118be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell" 24218be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 24318be62a5SSatish Balay { 24418be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 24518be62a5SSatish Balay PetscErrorCode ierr; 246*8fe8eb27SJed Brown Vec xx; 24718be62a5SSatish Balay 24818be62a5SSatish Balay PetscFunctionBegin; 249*8fe8eb27SJed Brown ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr); 250*8fe8eb27SJed Brown ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr); 251*8fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 252*8fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 25318be62a5SSatish Balay PetscFunctionReturn(0); 25418be62a5SSatish Balay } 25518be62a5SSatish Balay 25618be62a5SSatish Balay #undef __FUNCT__ 25718be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell" 25818be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 25918be62a5SSatish Balay { 26018be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 26118be62a5SSatish Balay PetscErrorCode ierr; 26218be62a5SSatish Balay 26318be62a5SSatish Balay PetscFunctionBegin; 26418be62a5SSatish Balay ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr); 26518be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 266*8fe8eb27SJed Brown if (shell->dshift) { 267*8fe8eb27SJed Brown ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr); 268*8fe8eb27SJed Brown } else { 26918be62a5SSatish Balay ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 27018be62a5SSatish Balay } 271*8fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 272*8fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 27318be62a5SSatish Balay PetscFunctionReturn(0); 27418be62a5SSatish Balay } 27518be62a5SSatish Balay 27618be62a5SSatish Balay #undef __FUNCT__ 277ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell" 278f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 279ef66eb69SBarry Smith { 280ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 281*8fe8eb27SJed Brown PetscErrorCode ierr; 282b24ad042SBarry Smith 283ef66eb69SBarry Smith PetscFunctionBegin; 284*8fe8eb27SJed Brown if (shell->left || shell->right || shell->dshift) { 285*8fe8eb27SJed Brown if (!shell->dshift) { 286*8fe8eb27SJed Brown if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);} 287*8fe8eb27SJed Brown shell->dshift = shell->dshift_owned; 288*8fe8eb27SJed Brown ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr); 289*8fe8eb27SJed Brown } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);} 290*8fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 291*8fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 292*8fe8eb27SJed Brown } else shell->vshift += a; 293*8fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 294ef66eb69SBarry Smith PetscFunctionReturn(0); 295ef66eb69SBarry Smith } 296ef66eb69SBarry Smith 297ef66eb69SBarry Smith #undef __FUNCT__ 298ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell" 299f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 300ef66eb69SBarry Smith { 301ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 302*8fe8eb27SJed Brown PetscErrorCode ierr; 303b24ad042SBarry Smith 304ef66eb69SBarry Smith PetscFunctionBegin; 305f4df32b1SMatthew Knepley shell->vscale *= a; 306*8fe8eb27SJed Brown if (shell->dshift) {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);} 307*8fe8eb27SJed Brown else shell->vshift *= a; 308*8fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 309*8fe8eb27SJed Brown PetscFunctionReturn(0); 31018be62a5SSatish Balay } 311*8fe8eb27SJed Brown 312*8fe8eb27SJed Brown #undef __FUNCT__ 313*8fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell" 314*8fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 315*8fe8eb27SJed Brown { 316*8fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 317*8fe8eb27SJed Brown PetscErrorCode ierr; 318*8fe8eb27SJed Brown 319*8fe8eb27SJed Brown PetscFunctionBegin; 320*8fe8eb27SJed Brown if (left) { 321*8fe8eb27SJed Brown if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);} 322*8fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);} 323*8fe8eb27SJed Brown else { 324*8fe8eb27SJed Brown shell->left = shell->left_owned; 325*8fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 32618be62a5SSatish Balay } 327ef66eb69SBarry Smith } 328*8fe8eb27SJed Brown if (right) { 329*8fe8eb27SJed Brown if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);} 330*8fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);} 331*8fe8eb27SJed Brown else { 332*8fe8eb27SJed Brown shell->right = shell->right_owned; 333*8fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 334*8fe8eb27SJed Brown } 335*8fe8eb27SJed Brown } 336*8fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 337ef66eb69SBarry Smith PetscFunctionReturn(0); 338ef66eb69SBarry Smith } 339ef66eb69SBarry Smith 340ef66eb69SBarry Smith #undef __FUNCT__ 341ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell" 342dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 343ef66eb69SBarry Smith { 344ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 345ef66eb69SBarry Smith 346ef66eb69SBarry Smith PetscFunctionBegin; 347*8fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 348ef66eb69SBarry Smith shell->vshift = 0.0; 349ef66eb69SBarry Smith shell->vscale = 1.0; 350*8fe8eb27SJed Brown shell->dshift = PETSC_NULL; 351*8fe8eb27SJed Brown shell->left = PETSC_NULL; 352*8fe8eb27SJed Brown shell->right = PETSC_NULL; 353*8fe8eb27SJed Brown shell->usingscaled = PETSC_FALSE; 354ef66eb69SBarry Smith Y->ops->mult = shell->mult; 35518be62a5SSatish Balay Y->ops->multtranspose = shell->multtranspose; 35618be62a5SSatish Balay Y->ops->getdiagonal = shell->getdiagonal; 357ef66eb69SBarry Smith } 358ef66eb69SBarry Smith PetscFunctionReturn(0); 359ef66eb69SBarry Smith } 360ef66eb69SBarry Smith 36109573ac7SBarry Smith extern PetscErrorCode MatConvert_Shell(Mat, const MatType,MatReuse,Mat*); 362b951964fSBarry Smith 363521d7252SBarry Smith #undef __FUNCT__ 364521d7252SBarry Smith #define __FUNCT__ "MatSetBlockSize_Shell" 365521d7252SBarry Smith PetscErrorCode MatSetBlockSize_Shell(Mat A,PetscInt bs) 366521d7252SBarry Smith { 36741c166b1SJed Brown PetscErrorCode ierr; 36841c166b1SJed Brown 369521d7252SBarry Smith PetscFunctionBegin; 37041c166b1SJed Brown ierr = PetscLayoutSetBlockSize(A->rmap,bs);CHKERRQ(ierr); 37141c166b1SJed Brown ierr = PetscLayoutSetBlockSize(A->cmap,bs);CHKERRQ(ierr); 372521d7252SBarry Smith PetscFunctionReturn(0); 373521d7252SBarry Smith } 374521d7252SBarry Smith 37509dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 37620563c6bSBarry Smith 0, 37720563c6bSBarry Smith 0, 37820563c6bSBarry Smith 0, 37997304618SKris Buschelman /* 4*/ 0, 38020563c6bSBarry Smith 0, 381b951964fSBarry Smith 0, 382b951964fSBarry Smith 0, 383b951964fSBarry Smith 0, 384b951964fSBarry Smith 0, 38597304618SKris Buschelman /*10*/ 0, 386b951964fSBarry Smith 0, 387b951964fSBarry Smith 0, 388b951964fSBarry Smith 0, 389b951964fSBarry Smith 0, 39097304618SKris Buschelman /*15*/ 0, 391b951964fSBarry Smith 0, 392b951964fSBarry Smith 0, 393*8fe8eb27SJed Brown MatDiagonalScale_Shell, 394b951964fSBarry Smith 0, 39597304618SKris Buschelman /*20*/ 0, 396ef66eb69SBarry Smith MatAssemblyEnd_Shell, 397b951964fSBarry Smith 0, 398b951964fSBarry Smith 0, 399d519adbfSMatthew Knepley /*24*/ 0, 400b951964fSBarry Smith 0, 401b951964fSBarry Smith 0, 402b951964fSBarry Smith 0, 403b951964fSBarry Smith 0, 404d519adbfSMatthew Knepley /*29*/ 0, 405b951964fSBarry Smith 0, 406273d9f13SBarry Smith 0, 407b951964fSBarry Smith 0, 408b951964fSBarry Smith 0, 409d519adbfSMatthew Knepley /*34*/ 0, 410b951964fSBarry Smith 0, 411b951964fSBarry Smith 0, 41209dc0095SBarry Smith 0, 41309dc0095SBarry Smith 0, 414d519adbfSMatthew Knepley /*39*/ 0, 41509dc0095SBarry Smith 0, 41609dc0095SBarry Smith 0, 41709dc0095SBarry Smith 0, 41809dc0095SBarry Smith 0, 419d519adbfSMatthew Knepley /*44*/ 0, 420ef66eb69SBarry Smith MatScale_Shell, 421ef66eb69SBarry Smith MatShift_Shell, 42209dc0095SBarry Smith 0, 42309dc0095SBarry Smith 0, 424d519adbfSMatthew Knepley /*49*/ MatSetBlockSize_Shell, 42509dc0095SBarry Smith 0, 42609dc0095SBarry Smith 0, 42709dc0095SBarry Smith 0, 42809dc0095SBarry Smith 0, 429d519adbfSMatthew Knepley /*54*/ 0, 43009dc0095SBarry Smith 0, 43109dc0095SBarry Smith 0, 43209dc0095SBarry Smith 0, 43309dc0095SBarry Smith 0, 434d519adbfSMatthew Knepley /*59*/ 0, 435b9b97703SBarry Smith MatDestroy_Shell, 43609dc0095SBarry Smith 0, 437357abbc8SBarry Smith 0, 438273d9f13SBarry Smith 0, 439d519adbfSMatthew Knepley /*64*/ 0, 440273d9f13SBarry Smith 0, 441273d9f13SBarry Smith 0, 442273d9f13SBarry Smith 0, 443273d9f13SBarry Smith 0, 444d519adbfSMatthew Knepley /*69*/ 0, 445273d9f13SBarry Smith 0, 446c87e5d42SMatthew Knepley MatConvert_Shell, 447273d9f13SBarry Smith 0, 44897304618SKris Buschelman 0, 449d519adbfSMatthew Knepley /*74*/ 0, 45097304618SKris Buschelman 0, 45197304618SKris Buschelman 0, 45297304618SKris Buschelman 0, 45397304618SKris Buschelman 0, 454d519adbfSMatthew Knepley /*79*/ 0, 45597304618SKris Buschelman 0, 45697304618SKris Buschelman 0, 45797304618SKris Buschelman 0, 458865e5f61SKris Buschelman 0, 459d519adbfSMatthew Knepley /*84*/ 0, 460865e5f61SKris Buschelman 0, 461865e5f61SKris Buschelman 0, 462865e5f61SKris Buschelman 0, 463865e5f61SKris Buschelman 0, 464d519adbfSMatthew Knepley /*89*/ 0, 465865e5f61SKris Buschelman 0, 466865e5f61SKris Buschelman 0, 467865e5f61SKris Buschelman 0, 468865e5f61SKris Buschelman 0, 469d519adbfSMatthew Knepley /*94*/ 0, 470865e5f61SKris Buschelman 0, 471865e5f61SKris Buschelman 0, 472865e5f61SKris Buschelman 0}; 473273d9f13SBarry Smith 4740bad9183SKris Buschelman /*MC 475fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 4760bad9183SKris Buschelman 4770bad9183SKris Buschelman Level: advanced 4780bad9183SKris Buschelman 4790bad9183SKris Buschelman .seealso: MatCreateShell 4800bad9183SKris Buschelman M*/ 4810bad9183SKris Buschelman 482273d9f13SBarry Smith EXTERN_C_BEGIN 483711e205bSSatish Balay #undef __FUNCT__ 484711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell" 4857087cfbeSBarry Smith PetscErrorCode MatCreate_Shell(Mat A) 486273d9f13SBarry Smith { 487273d9f13SBarry Smith Mat_Shell *b; 488dfbe8321SBarry Smith PetscErrorCode ierr; 489273d9f13SBarry Smith 490273d9f13SBarry Smith PetscFunctionBegin; 491273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 492273d9f13SBarry Smith 49338f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_Shell,&b);CHKERRQ(ierr); 494273d9f13SBarry Smith A->data = (void*)b; 495273d9f13SBarry Smith 49626283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 49726283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 49826283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 49926283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 500273d9f13SBarry Smith 501273d9f13SBarry Smith b->ctx = 0; 502ef66eb69SBarry Smith b->vshift = 0.0; 503ef66eb69SBarry Smith b->vscale = 1.0; 504ef66eb69SBarry Smith b->mult = 0; 50518be62a5SSatish Balay b->multtranspose = 0; 50618be62a5SSatish Balay b->getdiagonal = 0; 507273d9f13SBarry Smith A->assembled = PETSC_TRUE; 508210f0121SBarry Smith A->preallocated = PETSC_FALSE; 50917667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 510273d9f13SBarry Smith PetscFunctionReturn(0); 511273d9f13SBarry Smith } 512273d9f13SBarry Smith EXTERN_C_END 513e51e0e81SBarry Smith 514711e205bSSatish Balay #undef __FUNCT__ 515711e205bSSatish Balay #define __FUNCT__ "MatCreateShell" 5164b828684SBarry Smith /*@C 517052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 518ff756334SLois Curfman McInnes private data storage format. 519e51e0e81SBarry Smith 520c7fcc2eaSBarry Smith Collective on MPI_Comm 521c7fcc2eaSBarry Smith 522e51e0e81SBarry Smith Input Parameters: 523c7fcc2eaSBarry Smith + comm - MPI communicator 524c7fcc2eaSBarry Smith . m - number of local rows (must be given) 525c7fcc2eaSBarry Smith . n - number of local columns (must be given) 526c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 527c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 528c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 529e51e0e81SBarry Smith 530ff756334SLois Curfman McInnes Output Parameter: 53144cd7ae7SLois Curfman McInnes . A - the matrix 532e51e0e81SBarry Smith 533ff2fd236SBarry Smith Level: advanced 534ff2fd236SBarry Smith 535f39d1f56SLois Curfman McInnes Usage: 5367b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 537f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 538c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 539f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 540f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 541f39d1f56SLois Curfman McInnes 542ff756334SLois Curfman McInnes Notes: 543ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 544ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 545ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 546e51e0e81SBarry Smith 547f38a8d56SBarry Smith Fortran Notes: The context can only be an integer or a PetscObject 548f38a8d56SBarry Smith unfortunately it cannot be a Fortran array or derived type. 549f38a8d56SBarry Smith 550f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 551f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 552645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 553645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 554645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 555645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 556645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 557645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 558645985a0SLois Curfman McInnes For example, 559f39d1f56SLois Curfman McInnes 560f39d1f56SLois Curfman McInnes $ 561f39d1f56SLois Curfman McInnes $ Vec x, y 5627b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 563645985a0SLois Curfman McInnes $ Mat A 564f39d1f56SLois Curfman McInnes $ 565c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 566c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 567f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 568c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 569c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 570c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 571645985a0SLois Curfman McInnes $ MatMult(A,x,y); 572645985a0SLois Curfman McInnes $ MatDestroy(A); 573f39d1f56SLois Curfman McInnes $ VecDestroy(y); VecDestroy(x); 574645985a0SLois Curfman McInnes $ 575e51e0e81SBarry Smith 5760b627109SLois Curfman McInnes .keywords: matrix, shell, create 5770b627109SLois Curfman McInnes 578ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext() 579e51e0e81SBarry Smith @*/ 5807087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 581e51e0e81SBarry Smith { 582dfbe8321SBarry Smith PetscErrorCode ierr; 583ed3cc1f0SBarry Smith 5843a40ed3dSBarry Smith PetscFunctionBegin; 585f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 586f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 587899cda47SBarry Smith 588273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 589273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 590273d9f13SBarry Smith PetscFunctionReturn(0); 591c7fcc2eaSBarry Smith } 592c7fcc2eaSBarry Smith 593711e205bSSatish Balay #undef __FUNCT__ 594711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext" 595c6866cfdSSatish Balay /*@ 596273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 597c7fcc2eaSBarry Smith 5983f9fe445SBarry Smith Logically Collective on Mat 599c7fcc2eaSBarry Smith 600273d9f13SBarry Smith Input Parameters: 601273d9f13SBarry Smith + mat - the shell matrix 602273d9f13SBarry Smith - ctx - the context 603273d9f13SBarry Smith 604273d9f13SBarry Smith Level: advanced 605273d9f13SBarry Smith 606f38a8d56SBarry Smith Fortran Notes: The context can only be an integer or a PetscObject 607f38a8d56SBarry Smith unfortunately it cannot be a Fortran array or derived type. 608273d9f13SBarry Smith 609273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 6100bc0a719SBarry Smith @*/ 6117087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 612273d9f13SBarry Smith { 613273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 614dfbe8321SBarry Smith PetscErrorCode ierr; 615ace3abfcSBarry Smith PetscBool flg; 616273d9f13SBarry Smith 617273d9f13SBarry Smith PetscFunctionBegin; 6180700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 619273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 620273d9f13SBarry Smith if (flg) { 621273d9f13SBarry Smith shell->ctx = ctx; 622273d9f13SBarry Smith } 6233a40ed3dSBarry Smith PetscFunctionReturn(0); 624e51e0e81SBarry Smith } 625e51e0e81SBarry Smith 626711e205bSSatish Balay #undef __FUNCT__ 627711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation" 628c16cb8f2SBarry Smith /*@C 6293a3eedf2SBarry Smith MatShellSetOperation - Allows user to set a matrix operation for 6303a3eedf2SBarry Smith a shell matrix. 631e51e0e81SBarry Smith 6323f9fe445SBarry Smith Logically Collective on Mat 633fee21e36SBarry Smith 634c7fcc2eaSBarry Smith Input Parameters: 635c7fcc2eaSBarry Smith + mat - the shell matrix 636c7fcc2eaSBarry Smith . op - the name of the operation 637c7fcc2eaSBarry Smith - f - the function that provides the operation. 638c7fcc2eaSBarry Smith 63915091d37SBarry Smith Level: advanced 64015091d37SBarry Smith 641fae171e0SBarry Smith Usage: 642e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 643f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 644c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 6450b627109SLois Curfman McInnes 646a62d957aSLois Curfman McInnes Notes: 647e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 6481c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 649a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 6501c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 651a62d957aSLois Curfman McInnes 652a62d957aSLois Curfman McInnes All user-provided functions should have the same calling 653deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 654deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 655deebb3c3SLois Curfman McInnes routines, e.g., 656a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 657a62d957aSLois Curfman McInnes 6584aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 6594aa34b0aSBarry Smith nonzero on failure. 6604aa34b0aSBarry Smith 661a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 662a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 663a62d957aSLois Curfman McInnes set by MatCreateShell(). 664a62d957aSLois Curfman McInnes 665501d9185SBarry Smith Fortran Notes: For MatGetVecs() the user code should check if the input left or right matrix is -1 and in that case not 666501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 667501d9185SBarry Smith 668a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation 669a62d957aSLois Curfman McInnes 670ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext() 671e51e0e81SBarry Smith @*/ 6727087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 673e51e0e81SBarry Smith { 674dfbe8321SBarry Smith PetscErrorCode ierr; 675ace3abfcSBarry Smith PetscBool flg; 676273d9f13SBarry Smith 6773a40ed3dSBarry Smith PetscFunctionBegin; 6780700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6791c1c02c0SLois Curfman McInnes if (op == MATOP_DESTROY) { 680273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 681273d9f13SBarry Smith if (flg) { 682a62d957aSLois Curfman McInnes Mat_Shell *shell = (Mat_Shell*)mat->data; 6836849ba73SBarry Smith shell->destroy = (PetscErrorCode (*)(Mat)) f; 6846849ba73SBarry Smith } else mat->ops->destroy = (PetscErrorCode (*)(Mat)) f; 685a62d957aSLois Curfman McInnes } 6866849ba73SBarry Smith else if (op == MATOP_VIEW) mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer)) f; 687c134de8dSSatish Balay else (((void(**)(void))mat->ops)[op]) = f; 688a62d957aSLois Curfman McInnes 6893a40ed3dSBarry Smith PetscFunctionReturn(0); 690e51e0e81SBarry Smith } 691f0479e8cSBarry Smith 692711e205bSSatish Balay #undef __FUNCT__ 693711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation" 694d4bb536fSBarry Smith /*@C 695d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 696d4bb536fSBarry Smith 697c7fcc2eaSBarry Smith Not Collective 698c7fcc2eaSBarry Smith 699d4bb536fSBarry Smith Input Parameters: 700c7fcc2eaSBarry Smith + mat - the shell matrix 701c7fcc2eaSBarry Smith - op - the name of the operation 702d4bb536fSBarry Smith 703d4bb536fSBarry Smith Output Parameter: 704d4bb536fSBarry Smith . f - the function that provides the operation. 705d4bb536fSBarry Smith 70615091d37SBarry Smith Level: advanced 70715091d37SBarry Smith 708d4bb536fSBarry Smith Notes: 709e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 710d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 711d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 712d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 713d4bb536fSBarry Smith 714d4bb536fSBarry Smith All user-provided functions have the same calling 715d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 716d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 717d4bb536fSBarry Smith routines, e.g., 718d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 719d4bb536fSBarry Smith 720d4bb536fSBarry Smith Within each user-defined routine, the user should call 721d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 722d4bb536fSBarry Smith set by MatCreateShell(). 723d4bb536fSBarry Smith 724d4bb536fSBarry Smith .keywords: matrix, shell, set, operation 725d4bb536fSBarry Smith 726ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 727d4bb536fSBarry Smith @*/ 7287087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 729d4bb536fSBarry Smith { 730dfbe8321SBarry Smith PetscErrorCode ierr; 731ace3abfcSBarry Smith PetscBool flg; 732273d9f13SBarry Smith 7333a40ed3dSBarry Smith PetscFunctionBegin; 7340700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 735d4bb536fSBarry Smith if (op == MATOP_DESTROY) { 736273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 737273d9f13SBarry Smith if (flg) { 738d4bb536fSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 739c134de8dSSatish Balay *f = (void(*)(void))shell->destroy; 740c7fcc2eaSBarry Smith } else { 741c134de8dSSatish Balay *f = (void(*)(void))mat->ops->destroy; 742d4bb536fSBarry Smith } 743c7fcc2eaSBarry Smith } else if (op == MATOP_VIEW) { 744c134de8dSSatish Balay *f = (void(*)(void))mat->ops->view; 745c7fcc2eaSBarry Smith } else { 746c134de8dSSatish Balay *f = (((void(**)(void))mat->ops)[op]); 747d4bb536fSBarry Smith } 748d4bb536fSBarry Smith 7493a40ed3dSBarry Smith PetscFunctionReturn(0); 750d4bb536fSBarry Smith } 751d4bb536fSBarry Smith 752