1be1d678aSKris Buschelman 2e51e0e81SBarry Smith /* 320563c6bSBarry Smith This provides a simple shell for Fortran (and C programmers) to 420563c6bSBarry Smith create a very simple matrix class for use with KSP without coding 5ed3cc1f0SBarry Smith much of anything. 6e51e0e81SBarry Smith */ 7e51e0e81SBarry Smith 8af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 9e51e0e81SBarry Smith 1020563c6bSBarry Smith typedef struct { 116849ba73SBarry Smith PetscErrorCode (*destroy)(Mat); 126849ba73SBarry Smith PetscErrorCode (*mult)(Mat,Vec,Vec); 1318be62a5SSatish Balay PetscErrorCode (*multtranspose)(Mat,Vec,Vec); 1418be62a5SSatish Balay PetscErrorCode (*getdiagonal)(Mat,Vec); 15*6464bf51SAlex Fikl PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode); 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; 238e3f487b0SBarry Smith PetscObjectState instate,outstate; 239ef66eb69SBarry Smith 240ef66eb69SBarry Smith PetscFunctionBegin; 2418fe8eb27SJed Brown ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr); 242e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 2438fe8eb27SJed Brown ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr); 244e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 245e3f487b0SBarry Smith if (instate == outstate) { 246e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 247e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 248e3f487b0SBarry Smith } 2498fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2508fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 251ef66eb69SBarry Smith PetscFunctionReturn(0); 252ef66eb69SBarry Smith } 253ef66eb69SBarry Smith 254ef66eb69SBarry Smith #undef __FUNCT__ 2555edf6516SJed Brown #define __FUNCT__ "MatMultAdd_Shell" 2565edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 2575edf6516SJed Brown { 2585edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2595edf6516SJed Brown PetscErrorCode ierr; 2605edf6516SJed Brown 2615edf6516SJed Brown PetscFunctionBegin; 2625edf6516SJed Brown if (y == z) { 2635edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 2645edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 265b8430d49SBarry Smith ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr); 2665edf6516SJed Brown } else { 2675edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 2685edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 2695edf6516SJed Brown } 2705edf6516SJed Brown PetscFunctionReturn(0); 2715edf6516SJed Brown } 2725edf6516SJed Brown 2735edf6516SJed Brown #undef __FUNCT__ 27418be62a5SSatish Balay #define __FUNCT__ "MatMultTranspose_Shell" 27518be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 27618be62a5SSatish Balay { 27718be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 27818be62a5SSatish Balay PetscErrorCode ierr; 27925578ef6SJed Brown Vec xx; 280e3f487b0SBarry Smith PetscObjectState instate,outstate; 28118be62a5SSatish Balay 28218be62a5SSatish Balay PetscFunctionBegin; 2838fe8eb27SJed Brown ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr); 284e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 2858fe8eb27SJed Brown ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr); 286e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 287e3f487b0SBarry Smith if (instate == outstate) { 288e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 289e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 290e3f487b0SBarry Smith } 2918fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2928fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 29318be62a5SSatish Balay PetscFunctionReturn(0); 29418be62a5SSatish Balay } 29518be62a5SSatish Balay 29618be62a5SSatish Balay #undef __FUNCT__ 2975edf6516SJed Brown #define __FUNCT__ "MatMultTransposeAdd_Shell" 2985edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 2995edf6516SJed Brown { 3005edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 3015edf6516SJed Brown PetscErrorCode ierr; 3025edf6516SJed Brown 3035edf6516SJed Brown PetscFunctionBegin; 3045edf6516SJed Brown if (y == z) { 3055edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 3065edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 3075edf6516SJed Brown ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr); 3085edf6516SJed Brown } else { 3095edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 3105edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 3115edf6516SJed Brown } 3125edf6516SJed Brown PetscFunctionReturn(0); 3135edf6516SJed Brown } 3145edf6516SJed Brown 3155edf6516SJed Brown #undef __FUNCT__ 31618be62a5SSatish Balay #define __FUNCT__ "MatGetDiagonal_Shell" 31718be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 31818be62a5SSatish Balay { 31918be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 32018be62a5SSatish Balay PetscErrorCode ierr; 32118be62a5SSatish Balay 32218be62a5SSatish Balay PetscFunctionBegin; 32318be62a5SSatish Balay ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr); 32418be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 3258fe8eb27SJed Brown if (shell->dshift) { 3268fe8eb27SJed Brown ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr); 3278fe8eb27SJed Brown } else { 32818be62a5SSatish Balay ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 32918be62a5SSatish Balay } 3308fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 3318fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 33218be62a5SSatish Balay PetscFunctionReturn(0); 33318be62a5SSatish Balay } 33418be62a5SSatish Balay 33518be62a5SSatish Balay #undef __FUNCT__ 336ef66eb69SBarry Smith #define __FUNCT__ "MatShift_Shell" 337f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 338ef66eb69SBarry Smith { 339ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3408fe8eb27SJed Brown PetscErrorCode ierr; 341b24ad042SBarry Smith 342ef66eb69SBarry Smith PetscFunctionBegin; 3438fe8eb27SJed Brown if (shell->left || shell->right || shell->dshift) { 3448fe8eb27SJed Brown if (!shell->dshift) { 3458fe8eb27SJed Brown if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);} 3468fe8eb27SJed Brown shell->dshift = shell->dshift_owned; 3478fe8eb27SJed Brown ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr); 3488fe8eb27SJed Brown } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);} 3498fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 3508fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 3518fe8eb27SJed Brown } else shell->vshift += a; 3528fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 353ef66eb69SBarry Smith PetscFunctionReturn(0); 354ef66eb69SBarry Smith } 355ef66eb69SBarry Smith 356ef66eb69SBarry Smith #undef __FUNCT__ 357*6464bf51SAlex Fikl #define __FUNCT__ "MatDiagonalSet_Shell" 358*6464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 359*6464bf51SAlex Fikl { 360*6464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 361*6464bf51SAlex Fikl PetscErrorCode ierr; 362*6464bf51SAlex Fikl 363*6464bf51SAlex Fikl PetscFunctionBegin; 364*6464bf51SAlex Fikl if (shell->vscale != 1.0 || shell->left || shell->right) { 365*6464bf51SAlex Fikl SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported on scaled matrix"); 366*6464bf51SAlex Fikl } 367*6464bf51SAlex Fikl 368*6464bf51SAlex Fikl if (shell->diagonalset) { ierr = (*shell->diagonalset)(A,D,ins);CHKERRQ(ierr); } 369*6464bf51SAlex Fikl shell->vshift = 0.0; 370*6464bf51SAlex Fikl if (shell->dshift) { ierr = VecZeroEntries(shell->dshift);CHKERRQ(ierr); } 371*6464bf51SAlex Fikl PetscFunctionReturn(0); 372*6464bf51SAlex Fikl } 373*6464bf51SAlex Fikl 374*6464bf51SAlex Fikl #undef __FUNCT__ 375ef66eb69SBarry Smith #define __FUNCT__ "MatScale_Shell" 376f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 377ef66eb69SBarry Smith { 378ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3798fe8eb27SJed Brown PetscErrorCode ierr; 380b24ad042SBarry Smith 381ef66eb69SBarry Smith PetscFunctionBegin; 382f4df32b1SMatthew Knepley shell->vscale *= a; 3832205254eSKarl Rupp if (shell->dshift) { 3842205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 3852205254eSKarl Rupp } else shell->vshift *= a; 3868fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 3878fe8eb27SJed Brown PetscFunctionReturn(0); 38818be62a5SSatish Balay } 3898fe8eb27SJed Brown 3908fe8eb27SJed Brown #undef __FUNCT__ 3918fe8eb27SJed Brown #define __FUNCT__ "MatDiagonalScale_Shell" 3928fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 3938fe8eb27SJed Brown { 3948fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 3958fe8eb27SJed Brown PetscErrorCode ierr; 3968fe8eb27SJed Brown 3978fe8eb27SJed Brown PetscFunctionBegin; 3988fe8eb27SJed Brown if (left) { 3998fe8eb27SJed Brown if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);} 4002205254eSKarl Rupp if (shell->left) { 4012205254eSKarl Rupp ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 4022205254eSKarl Rupp } else { 4038fe8eb27SJed Brown shell->left = shell->left_owned; 4048fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 40518be62a5SSatish Balay } 406ef66eb69SBarry Smith } 4078fe8eb27SJed Brown if (right) { 4088fe8eb27SJed Brown if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);} 4092205254eSKarl Rupp if (shell->right) { 4102205254eSKarl Rupp ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 4112205254eSKarl Rupp } else { 4128fe8eb27SJed Brown shell->right = shell->right_owned; 4138fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 4148fe8eb27SJed Brown } 4158fe8eb27SJed Brown } 4168fe8eb27SJed Brown ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr); 417ef66eb69SBarry Smith PetscFunctionReturn(0); 418ef66eb69SBarry Smith } 419ef66eb69SBarry Smith 420ef66eb69SBarry Smith #undef __FUNCT__ 421ef66eb69SBarry Smith #define __FUNCT__ "MatAssemblyEnd_Shell" 422dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 423ef66eb69SBarry Smith { 424ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 425ef66eb69SBarry Smith 426ef66eb69SBarry Smith PetscFunctionBegin; 4278fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 428ef66eb69SBarry Smith shell->vshift = 0.0; 429ef66eb69SBarry Smith shell->vscale = 1.0; 4300298fd71SBarry Smith shell->dshift = NULL; 4310298fd71SBarry Smith shell->left = NULL; 4320298fd71SBarry Smith shell->right = NULL; 4337fb7d96aSJed Brown if (shell->mult) { 434ef66eb69SBarry Smith Y->ops->mult = shell->mult; 4350298fd71SBarry Smith shell->mult = NULL; 4367fb7d96aSJed Brown } 4377fb7d96aSJed Brown if (shell->multtranspose) { 43818be62a5SSatish Balay Y->ops->multtranspose = shell->multtranspose; 4390298fd71SBarry Smith shell->multtranspose = NULL; 4407fb7d96aSJed Brown } 4417fb7d96aSJed Brown if (shell->getdiagonal) { 44218be62a5SSatish Balay Y->ops->getdiagonal = shell->getdiagonal; 4430298fd71SBarry Smith shell->getdiagonal = NULL; 4447fb7d96aSJed Brown } 4457fb7d96aSJed Brown shell->usingscaled = PETSC_FALSE; 446ef66eb69SBarry Smith } 447ef66eb69SBarry Smith PetscFunctionReturn(0); 448ef66eb69SBarry Smith } 449ef66eb69SBarry Smith 450cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*); 451b951964fSBarry Smith 4523b49f96aSBarry Smith #undef __FUNCT__ 4533b49f96aSBarry Smith #define __FUNCT__ "MatMissingDiagonal_Shell" 4543b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 4553b49f96aSBarry Smith { 4563b49f96aSBarry Smith PetscFunctionBegin; 4573b49f96aSBarry Smith *missing = PETSC_FALSE; 4583b49f96aSBarry Smith PetscFunctionReturn(0); 4593b49f96aSBarry Smith } 4603b49f96aSBarry Smith 46109dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 46220563c6bSBarry Smith 0, 46320563c6bSBarry Smith 0, 46420563c6bSBarry Smith 0, 46597304618SKris Buschelman /* 4*/ 0, 46620563c6bSBarry Smith 0, 467b951964fSBarry Smith 0, 468b951964fSBarry Smith 0, 469b951964fSBarry Smith 0, 470b951964fSBarry Smith 0, 47197304618SKris Buschelman /*10*/ 0, 472b951964fSBarry Smith 0, 473b951964fSBarry Smith 0, 474b951964fSBarry Smith 0, 475b951964fSBarry Smith 0, 47697304618SKris Buschelman /*15*/ 0, 477b951964fSBarry Smith 0, 478b951964fSBarry Smith 0, 4798fe8eb27SJed Brown MatDiagonalScale_Shell, 480b951964fSBarry Smith 0, 48197304618SKris Buschelman /*20*/ 0, 482ef66eb69SBarry Smith MatAssemblyEnd_Shell, 483b951964fSBarry Smith 0, 484b951964fSBarry Smith 0, 485d519adbfSMatthew Knepley /*24*/ 0, 486b951964fSBarry Smith 0, 487b951964fSBarry Smith 0, 488b951964fSBarry Smith 0, 489b951964fSBarry Smith 0, 490d519adbfSMatthew Knepley /*29*/ 0, 491b951964fSBarry Smith 0, 492273d9f13SBarry Smith 0, 493b951964fSBarry Smith 0, 494b951964fSBarry Smith 0, 495d519adbfSMatthew Knepley /*34*/ 0, 496b951964fSBarry Smith 0, 497b951964fSBarry Smith 0, 49809dc0095SBarry Smith 0, 49909dc0095SBarry Smith 0, 500d519adbfSMatthew Knepley /*39*/ 0, 50109dc0095SBarry Smith 0, 50209dc0095SBarry Smith 0, 50309dc0095SBarry Smith 0, 50409dc0095SBarry Smith 0, 505d519adbfSMatthew Knepley /*44*/ 0, 506ef66eb69SBarry Smith MatScale_Shell, 507ef66eb69SBarry Smith MatShift_Shell, 508*6464bf51SAlex Fikl MatDiagonalSet_Shell, 50909dc0095SBarry Smith 0, 510f73d5cc4SBarry Smith /*49*/ 0, 51109dc0095SBarry Smith 0, 51209dc0095SBarry Smith 0, 51309dc0095SBarry Smith 0, 51409dc0095SBarry Smith 0, 515d519adbfSMatthew Knepley /*54*/ 0, 51609dc0095SBarry Smith 0, 51709dc0095SBarry Smith 0, 51809dc0095SBarry Smith 0, 51909dc0095SBarry Smith 0, 520d519adbfSMatthew Knepley /*59*/ 0, 521b9b97703SBarry Smith MatDestroy_Shell, 52209dc0095SBarry Smith 0, 523357abbc8SBarry Smith 0, 524273d9f13SBarry Smith 0, 525d519adbfSMatthew Knepley /*64*/ 0, 526273d9f13SBarry Smith 0, 527273d9f13SBarry Smith 0, 528273d9f13SBarry Smith 0, 529273d9f13SBarry Smith 0, 530d519adbfSMatthew Knepley /*69*/ 0, 531273d9f13SBarry Smith 0, 532c87e5d42SMatthew Knepley MatConvert_Shell, 533273d9f13SBarry Smith 0, 53497304618SKris Buschelman 0, 535d519adbfSMatthew Knepley /*74*/ 0, 53697304618SKris Buschelman 0, 53797304618SKris Buschelman 0, 53897304618SKris Buschelman 0, 53997304618SKris Buschelman 0, 540d519adbfSMatthew Knepley /*79*/ 0, 54197304618SKris Buschelman 0, 54297304618SKris Buschelman 0, 54397304618SKris Buschelman 0, 544865e5f61SKris Buschelman 0, 545d519adbfSMatthew Knepley /*84*/ 0, 546865e5f61SKris Buschelman 0, 547865e5f61SKris Buschelman 0, 548865e5f61SKris Buschelman 0, 549865e5f61SKris Buschelman 0, 550d519adbfSMatthew Knepley /*89*/ 0, 551865e5f61SKris Buschelman 0, 552865e5f61SKris Buschelman 0, 553865e5f61SKris Buschelman 0, 554865e5f61SKris Buschelman 0, 555d519adbfSMatthew Knepley /*94*/ 0, 556865e5f61SKris Buschelman 0, 557865e5f61SKris Buschelman 0, 5583964eb88SJed Brown 0, 5593964eb88SJed Brown 0, 5603964eb88SJed Brown /*99*/ 0, 5613964eb88SJed Brown 0, 5623964eb88SJed Brown 0, 5633964eb88SJed Brown 0, 5643964eb88SJed Brown 0, 5653964eb88SJed Brown /*104*/ 0, 5663964eb88SJed Brown 0, 5673964eb88SJed Brown 0, 5683964eb88SJed Brown 0, 5693964eb88SJed Brown 0, 5703964eb88SJed Brown /*109*/ 0, 5713964eb88SJed Brown 0, 5723964eb88SJed Brown 0, 5733964eb88SJed Brown 0, 5743b49f96aSBarry Smith MatMissingDiagonal_Shell, 5753964eb88SJed Brown /*114*/ 0, 5763964eb88SJed Brown 0, 5773964eb88SJed Brown 0, 5783964eb88SJed Brown 0, 5793964eb88SJed Brown 0, 5803964eb88SJed Brown /*119*/ 0, 5813964eb88SJed Brown 0, 5823964eb88SJed Brown 0, 5833964eb88SJed Brown 0, 5843964eb88SJed Brown 0, 5853964eb88SJed Brown /*124*/ 0, 5863964eb88SJed Brown 0, 5873964eb88SJed Brown 0, 5883964eb88SJed Brown 0, 5893964eb88SJed Brown 0, 5903964eb88SJed Brown /*129*/ 0, 5913964eb88SJed Brown 0, 5923964eb88SJed Brown 0, 5933964eb88SJed Brown 0, 5943964eb88SJed Brown 0, 5953964eb88SJed Brown /*134*/ 0, 5963964eb88SJed Brown 0, 5973964eb88SJed Brown 0, 5983964eb88SJed Brown 0, 5993964eb88SJed Brown 0, 6003964eb88SJed Brown /*139*/ 0, 601f9426fe0SMark Adams 0, 6023964eb88SJed Brown 0 6033964eb88SJed Brown }; 604273d9f13SBarry Smith 6050bad9183SKris Buschelman /*MC 606fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 6070bad9183SKris Buschelman 6080bad9183SKris Buschelman Level: advanced 6090bad9183SKris Buschelman 6100bad9183SKris Buschelman .seealso: MatCreateShell 6110bad9183SKris Buschelman M*/ 6120bad9183SKris Buschelman 613711e205bSSatish Balay #undef __FUNCT__ 614711e205bSSatish Balay #define __FUNCT__ "MatCreate_Shell" 6158cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 616273d9f13SBarry Smith { 617273d9f13SBarry Smith Mat_Shell *b; 618dfbe8321SBarry Smith PetscErrorCode ierr; 619273d9f13SBarry Smith 620273d9f13SBarry Smith PetscFunctionBegin; 621273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 622273d9f13SBarry Smith 623b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 624273d9f13SBarry Smith A->data = (void*)b; 625273d9f13SBarry Smith 62626283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 62726283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 628273d9f13SBarry Smith 629273d9f13SBarry Smith b->ctx = 0; 630ef66eb69SBarry Smith b->vshift = 0.0; 631ef66eb69SBarry Smith b->vscale = 1.0; 632ef66eb69SBarry Smith b->mult = 0; 63318be62a5SSatish Balay b->multtranspose = 0; 63418be62a5SSatish Balay b->getdiagonal = 0; 635273d9f13SBarry Smith A->assembled = PETSC_TRUE; 636210f0121SBarry Smith A->preallocated = PETSC_FALSE; 6372205254eSKarl Rupp 63817667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 639273d9f13SBarry Smith PetscFunctionReturn(0); 640273d9f13SBarry Smith } 641e51e0e81SBarry Smith 642711e205bSSatish Balay #undef __FUNCT__ 643711e205bSSatish Balay #define __FUNCT__ "MatCreateShell" 6444b828684SBarry Smith /*@C 645052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 646ff756334SLois Curfman McInnes private data storage format. 647e51e0e81SBarry Smith 648c7fcc2eaSBarry Smith Collective on MPI_Comm 649c7fcc2eaSBarry Smith 650e51e0e81SBarry Smith Input Parameters: 651c7fcc2eaSBarry Smith + comm - MPI communicator 652c7fcc2eaSBarry Smith . m - number of local rows (must be given) 653c7fcc2eaSBarry Smith . n - number of local columns (must be given) 654c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 655c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 656c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 657e51e0e81SBarry Smith 658ff756334SLois Curfman McInnes Output Parameter: 65944cd7ae7SLois Curfman McInnes . A - the matrix 660e51e0e81SBarry Smith 661ff2fd236SBarry Smith Level: advanced 662ff2fd236SBarry Smith 663f39d1f56SLois Curfman McInnes Usage: 6647b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 665f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 666c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 667f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 668f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 669f39d1f56SLois Curfman McInnes 670ff756334SLois Curfman McInnes Notes: 671ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 672ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 673ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 674e51e0e81SBarry Smith 675daf670e6SBarry Smith Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this 676daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 677daf670e6SBarry Smith in as the ctx argument. 678f38a8d56SBarry Smith 679f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 680f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 681645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 682645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 683645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 684645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 685645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 686645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 687645985a0SLois Curfman McInnes For example, 688f39d1f56SLois Curfman McInnes 689f39d1f56SLois Curfman McInnes $ 690f39d1f56SLois Curfman McInnes $ Vec x, y 6917b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 692645985a0SLois Curfman McInnes $ Mat A 693f39d1f56SLois Curfman McInnes $ 694c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 695c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 696f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 697c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 698c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 699c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 700645985a0SLois Curfman McInnes $ MatMult(A,x,y); 701645985a0SLois Curfman McInnes $ MatDestroy(A); 702f39d1f56SLois Curfman McInnes $ VecDestroy(y); VecDestroy(x); 703645985a0SLois Curfman McInnes $ 704e51e0e81SBarry Smith 7050b627109SLois Curfman McInnes .keywords: matrix, shell, create 7060b627109SLois Curfman McInnes 707ab50ec6bSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext() 708e51e0e81SBarry Smith @*/ 7097087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 710e51e0e81SBarry Smith { 711dfbe8321SBarry Smith PetscErrorCode ierr; 712ed3cc1f0SBarry Smith 7133a40ed3dSBarry Smith PetscFunctionBegin; 714f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 715f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 716273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 717273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 7180fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 719273d9f13SBarry Smith PetscFunctionReturn(0); 720c7fcc2eaSBarry Smith } 721c7fcc2eaSBarry Smith 722711e205bSSatish Balay #undef __FUNCT__ 723711e205bSSatish Balay #define __FUNCT__ "MatShellSetContext" 724c6866cfdSSatish Balay /*@ 725273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 726c7fcc2eaSBarry Smith 7273f9fe445SBarry Smith Logically Collective on Mat 728c7fcc2eaSBarry Smith 729273d9f13SBarry Smith Input Parameters: 730273d9f13SBarry Smith + mat - the shell matrix 731273d9f13SBarry Smith - ctx - the context 732273d9f13SBarry Smith 733273d9f13SBarry Smith Level: advanced 734273d9f13SBarry Smith 735daf670e6SBarry Smith Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this 736daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 737273d9f13SBarry Smith 738273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 7390bc0a719SBarry Smith @*/ 7407087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 741273d9f13SBarry Smith { 742273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 743dfbe8321SBarry Smith PetscErrorCode ierr; 744ace3abfcSBarry Smith PetscBool flg; 745273d9f13SBarry Smith 746273d9f13SBarry Smith PetscFunctionBegin; 7470700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 748251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 749273d9f13SBarry Smith if (flg) { 750273d9f13SBarry Smith shell->ctx = ctx; 751ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix"); 7523a40ed3dSBarry Smith PetscFunctionReturn(0); 753e51e0e81SBarry Smith } 754e51e0e81SBarry Smith 755711e205bSSatish Balay #undef __FUNCT__ 756711e205bSSatish Balay #define __FUNCT__ "MatShellSetOperation" 757c16cb8f2SBarry Smith /*@C 7583a3eedf2SBarry Smith MatShellSetOperation - Allows user to set a matrix operation for 7593a3eedf2SBarry Smith a shell matrix. 760e51e0e81SBarry Smith 7613f9fe445SBarry Smith Logically Collective on Mat 762fee21e36SBarry Smith 763c7fcc2eaSBarry Smith Input Parameters: 764c7fcc2eaSBarry Smith + mat - the shell matrix 765c7fcc2eaSBarry Smith . op - the name of the operation 766c7fcc2eaSBarry Smith - f - the function that provides the operation. 767c7fcc2eaSBarry Smith 76815091d37SBarry Smith Level: advanced 76915091d37SBarry Smith 770fae171e0SBarry Smith Usage: 771e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 772f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 773c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 7740b627109SLois Curfman McInnes 775a62d957aSLois Curfman McInnes Notes: 776e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 7771c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 778a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 7791c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 780a62d957aSLois Curfman McInnes 78125296bd5SBarry Smith All user-provided functions (execept for MATOP_DESTROY) should have the same calling 782deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 783deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 784deebb3c3SLois Curfman McInnes routines, e.g., 785a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 786a62d957aSLois Curfman McInnes 7874aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 7884aa34b0aSBarry Smith nonzero on failure. 7894aa34b0aSBarry Smith 790a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 791a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 792a62d957aSLois Curfman McInnes set by MatCreateShell(). 793a62d957aSLois Curfman McInnes 7942a7a6963SBarry Smith Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 795501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 796501d9185SBarry Smith 797a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation 798a62d957aSLois Curfman McInnes 799ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext() 800e51e0e81SBarry Smith @*/ 8017087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 802e51e0e81SBarry Smith { 803dfbe8321SBarry Smith PetscErrorCode ierr; 804ace3abfcSBarry Smith PetscBool flg; 805273d9f13SBarry Smith 8063a40ed3dSBarry Smith PetscFunctionBegin; 8070700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 8085edf6516SJed Brown switch (op) { 8095edf6516SJed Brown case MATOP_DESTROY: 810251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 811273d9f13SBarry Smith if (flg) { 812a62d957aSLois Curfman McInnes Mat_Shell *shell = (Mat_Shell*)mat->data; 8136849ba73SBarry Smith shell->destroy = (PetscErrorCode (*)(Mat))f; 8146849ba73SBarry Smith } else mat->ops->destroy = (PetscErrorCode (*)(Mat))f; 8155edf6516SJed Brown break; 816*6464bf51SAlex Fikl case MATOP_DIAGONAL_SET: 817*6464bf51SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 818*6464bf51SAlex Fikl if (flg) { 819*6464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)mat->data; 820*6464bf51SAlex Fikl shell->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f; 821*6464bf51SAlex Fikl } else { 822*6464bf51SAlex Fikl mat->ops->diagonalset = (PetscErrorCode (*)(Mat,Vec,InsertMode))f; 823*6464bf51SAlex Fikl } 824*6464bf51SAlex Fikl break; 8255edf6516SJed Brown case MATOP_VIEW: 8265edf6516SJed Brown mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 8275edf6516SJed Brown break; 8285edf6516SJed Brown case MATOP_MULT: 8295edf6516SJed Brown mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 8305edf6516SJed Brown if (!mat->ops->multadd) mat->ops->multadd = MatMultAdd_Shell; 8315edf6516SJed Brown break; 8325edf6516SJed Brown case MATOP_MULT_TRANSPOSE: 8335edf6516SJed Brown mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 8345edf6516SJed Brown if (!mat->ops->multtransposeadd) mat->ops->multtransposeadd = MatMultTransposeAdd_Shell; 8355edf6516SJed Brown break; 8365edf6516SJed Brown default: 8375edf6516SJed Brown (((void(**)(void))mat->ops)[op]) = f; 838a62d957aSLois Curfman McInnes } 8393a40ed3dSBarry Smith PetscFunctionReturn(0); 840e51e0e81SBarry Smith } 841f0479e8cSBarry Smith 842711e205bSSatish Balay #undef __FUNCT__ 843711e205bSSatish Balay #define __FUNCT__ "MatShellGetOperation" 844d4bb536fSBarry Smith /*@C 845d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 846d4bb536fSBarry Smith 847c7fcc2eaSBarry Smith Not Collective 848c7fcc2eaSBarry Smith 849d4bb536fSBarry Smith Input Parameters: 850c7fcc2eaSBarry Smith + mat - the shell matrix 851c7fcc2eaSBarry Smith - op - the name of the operation 852d4bb536fSBarry Smith 853d4bb536fSBarry Smith Output Parameter: 854d4bb536fSBarry Smith . f - the function that provides the operation. 855d4bb536fSBarry Smith 85615091d37SBarry Smith Level: advanced 85715091d37SBarry Smith 858d4bb536fSBarry Smith Notes: 859e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 860d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 861d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 862d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 863d4bb536fSBarry Smith 864d4bb536fSBarry Smith All user-provided functions have the same calling 865d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 866d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 867d4bb536fSBarry Smith routines, e.g., 868d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 869d4bb536fSBarry Smith 870d4bb536fSBarry Smith Within each user-defined routine, the user should call 871d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 872d4bb536fSBarry Smith set by MatCreateShell(). 873d4bb536fSBarry Smith 874d4bb536fSBarry Smith .keywords: matrix, shell, set, operation 875d4bb536fSBarry Smith 876ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 877d4bb536fSBarry Smith @*/ 8787087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 879d4bb536fSBarry Smith { 880dfbe8321SBarry Smith PetscErrorCode ierr; 881ace3abfcSBarry Smith PetscBool flg; 882273d9f13SBarry Smith 8833a40ed3dSBarry Smith PetscFunctionBegin; 8840700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 885d4bb536fSBarry Smith if (op == MATOP_DESTROY) { 886251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 887273d9f13SBarry Smith if (flg) { 888d4bb536fSBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 889c134de8dSSatish Balay *f = (void (*)(void))shell->destroy; 890c7fcc2eaSBarry Smith } else { 891c134de8dSSatish Balay *f = (void (*)(void))mat->ops->destroy; 892d4bb536fSBarry Smith } 893c7fcc2eaSBarry Smith } else if (op == MATOP_VIEW) { 894c134de8dSSatish Balay *f = (void (*)(void))mat->ops->view; 895c7fcc2eaSBarry Smith } else { 896c134de8dSSatish Balay *f = (((void (**)(void))mat->ops)[op]); 897d4bb536fSBarry Smith } 8983a40ed3dSBarry Smith PetscFunctionReturn(0); 899d4bb536fSBarry Smith } 900d4bb536fSBarry Smith 901