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 1028f357e3SAlex Fikl struct _MatShellOps { 11976e8c5aSLisandro Dalcin /* 3 */ PetscErrorCode (*mult)(Mat,Vec,Vec); 12976e8c5aSLisandro Dalcin /* 5 */ PetscErrorCode (*multtranspose)(Mat,Vec,Vec); 13976e8c5aSLisandro Dalcin /* 17 */ PetscErrorCode (*getdiagonal)(Mat,Vec); 14976e8c5aSLisandro Dalcin /* 43 */ PetscErrorCode (*copy)(Mat,Mat,MatStructure); 15976e8c5aSLisandro Dalcin /* 60 */ PetscErrorCode (*destroy)(Mat); 1628f357e3SAlex Fikl }; 1728f357e3SAlex Fikl 1828f357e3SAlex Fikl typedef struct { 1928f357e3SAlex Fikl struct _MatShellOps ops[1]; 202205254eSKarl Rupp 21ef66eb69SBarry Smith PetscScalar vscale,vshift; 228fe8eb27SJed Brown Vec dshift; 238fe8eb27SJed Brown Vec left,right; 248fe8eb27SJed Brown Vec left_work,right_work; 255edf6516SJed Brown Vec left_add_work,right_add_work; 269f137db4SBarry Smith Mat axpy; 279f137db4SBarry Smith PetscScalar axpy_vscale; 280c0fd78eSBarry Smith PetscBool managescalingshifts; /* The user will manage the scaling and shifts for the MATSHELL, not the default */ 2920563c6bSBarry Smith void *ctx; 3088cf3e7dSBarry Smith } Mat_Shell; 310c0fd78eSBarry Smith 328fe8eb27SJed Brown /* 330c0fd78eSBarry Smith xx = diag(left)*x 348fe8eb27SJed Brown */ 358fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx) 368fe8eb27SJed Brown { 378fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 388fe8eb27SJed Brown PetscErrorCode ierr; 398fe8eb27SJed Brown 408fe8eb27SJed Brown PetscFunctionBegin; 410298fd71SBarry Smith *xx = NULL; 428fe8eb27SJed Brown if (!shell->left) { 438fe8eb27SJed Brown *xx = x; 448fe8eb27SJed Brown } else { 458fe8eb27SJed Brown if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);} 468fe8eb27SJed Brown ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr); 478fe8eb27SJed Brown *xx = shell->left_work; 488fe8eb27SJed Brown } 498fe8eb27SJed Brown PetscFunctionReturn(0); 508fe8eb27SJed Brown } 518fe8eb27SJed Brown 520c0fd78eSBarry Smith /* 530c0fd78eSBarry Smith xx = diag(right)*x 540c0fd78eSBarry Smith */ 558fe8eb27SJed Brown static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx) 568fe8eb27SJed Brown { 578fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 588fe8eb27SJed Brown PetscErrorCode ierr; 598fe8eb27SJed Brown 608fe8eb27SJed Brown PetscFunctionBegin; 610298fd71SBarry Smith *xx = NULL; 628fe8eb27SJed Brown if (!shell->right) { 638fe8eb27SJed Brown *xx = x; 648fe8eb27SJed Brown } else { 658fe8eb27SJed Brown if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);} 668fe8eb27SJed Brown ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr); 678fe8eb27SJed Brown *xx = shell->right_work; 688fe8eb27SJed Brown } 698fe8eb27SJed Brown PetscFunctionReturn(0); 708fe8eb27SJed Brown } 718fe8eb27SJed Brown 720c0fd78eSBarry Smith /* 730c0fd78eSBarry Smith x = diag(left)*x 740c0fd78eSBarry Smith */ 758fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x) 768fe8eb27SJed Brown { 778fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 788fe8eb27SJed Brown PetscErrorCode ierr; 798fe8eb27SJed Brown 808fe8eb27SJed Brown PetscFunctionBegin; 818fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);} 828fe8eb27SJed Brown PetscFunctionReturn(0); 838fe8eb27SJed Brown } 848fe8eb27SJed Brown 850c0fd78eSBarry Smith /* 860c0fd78eSBarry Smith x = diag(right)*x 870c0fd78eSBarry Smith */ 888fe8eb27SJed Brown static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x) 898fe8eb27SJed Brown { 908fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 918fe8eb27SJed Brown PetscErrorCode ierr; 928fe8eb27SJed Brown 938fe8eb27SJed Brown PetscFunctionBegin; 948fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);} 958fe8eb27SJed Brown PetscFunctionReturn(0); 968fe8eb27SJed Brown } 978fe8eb27SJed Brown 980c0fd78eSBarry Smith /* 990c0fd78eSBarry Smith Y = vscale*Y + diag(dshift)*X + vshift*X 1000c0fd78eSBarry Smith 1010c0fd78eSBarry Smith On input Y already contains A*x 1020c0fd78eSBarry Smith */ 1038fe8eb27SJed Brown static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y) 1048fe8eb27SJed Brown { 1058fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 1068fe8eb27SJed Brown PetscErrorCode ierr; 1078fe8eb27SJed Brown 1088fe8eb27SJed Brown PetscFunctionBegin; 1098fe8eb27SJed Brown if (shell->dshift) { /* get arrays because there is no VecPointwiseMultAdd() */ 1108fe8eb27SJed Brown PetscInt i,m; 1118fe8eb27SJed Brown const PetscScalar *x,*d; 1128fe8eb27SJed Brown PetscScalar *y; 1138fe8eb27SJed Brown ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr); 1148fe8eb27SJed Brown ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1158fe8eb27SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 1168fe8eb27SJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1178fe8eb27SJed Brown for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i]; 1188fe8eb27SJed Brown ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr); 1198fe8eb27SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 1208fe8eb27SJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 1210c0fd78eSBarry Smith } else { 122027c5fb5SJed Brown ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr); 1238fe8eb27SJed Brown } 124d4c7de66SBarry Smith if (shell->vshift != 0.0) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */ 1258fe8eb27SJed Brown PetscFunctionReturn(0); 1268fe8eb27SJed Brown } 127e51e0e81SBarry Smith 1289d225801SJed Brown /*@ 129a62d957aSLois Curfman McInnes MatShellGetContext - Returns the user-provided context associated with a shell matrix. 130b4fd4287SBarry Smith 131c7fcc2eaSBarry Smith Not Collective 132c7fcc2eaSBarry Smith 133b4fd4287SBarry Smith Input Parameter: 134b4fd4287SBarry Smith . mat - the matrix, should have been created with MatCreateShell() 135b4fd4287SBarry Smith 136b4fd4287SBarry Smith Output Parameter: 137b4fd4287SBarry Smith . ctx - the user provided context 138b4fd4287SBarry Smith 13915091d37SBarry Smith Level: advanced 14015091d37SBarry Smith 141*95452b02SPatrick Sanan Fortran Notes: 142*95452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 143daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 144a62d957aSLois Curfman McInnes 145a62d957aSLois Curfman McInnes .keywords: matrix, shell, get, context 146a62d957aSLois Curfman McInnes 147ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext() 1480bc0a719SBarry Smith @*/ 1498fe8eb27SJed Brown PetscErrorCode MatShellGetContext(Mat mat,void *ctx) 150b4fd4287SBarry Smith { 151dfbe8321SBarry Smith PetscErrorCode ierr; 152ace3abfcSBarry Smith PetscBool flg; 153273d9f13SBarry Smith 1543a40ed3dSBarry Smith PetscFunctionBegin; 1550700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1564482741eSBarry Smith PetscValidPointer(ctx,2); 157251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 158940183f3SBarry Smith if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx; 159ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix"); 1603a40ed3dSBarry Smith PetscFunctionReturn(0); 161b4fd4287SBarry Smith } 162b4fd4287SBarry Smith 163dfbe8321SBarry Smith PetscErrorCode MatDestroy_Shell(Mat mat) 164e51e0e81SBarry Smith { 165dfbe8321SBarry Smith PetscErrorCode ierr; 166bf0cc555SLisandro Dalcin Mat_Shell *shell = (Mat_Shell*)mat->data; 167ed3cc1f0SBarry Smith 1683a40ed3dSBarry Smith PetscFunctionBegin; 16928f357e3SAlex Fikl if (shell->ops->destroy) { 17028f357e3SAlex Fikl ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr); 171bf0cc555SLisandro Dalcin } 1720c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 1730c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 1740c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 1758fe8eb27SJed Brown ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr); 1768fe8eb27SJed Brown ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr); 1775edf6516SJed Brown ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr); 1785edf6516SJed Brown ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr); 1799f137db4SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 180bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 1813a40ed3dSBarry Smith PetscFunctionReturn(0); 182e51e0e81SBarry Smith } 183e51e0e81SBarry Smith 1847fabda3fSAlex Fikl PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str) 1857fabda3fSAlex Fikl { 18628f357e3SAlex Fikl Mat_Shell *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data; 1877fabda3fSAlex Fikl PetscErrorCode ierr; 188cb8c8a77SBarry Smith PetscBool matflg; 1897fabda3fSAlex Fikl 1907fabda3fSAlex Fikl PetscFunctionBegin; 19128f357e3SAlex Fikl ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr); 192cb8c8a77SBarry Smith if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL"); 1937fabda3fSAlex Fikl 194cb8c8a77SBarry Smith ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 195cb8c8a77SBarry Smith ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr); 1967fabda3fSAlex Fikl 197cb8c8a77SBarry Smith if (shellA->ops->copy) { 19828f357e3SAlex Fikl ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr); 199cb8c8a77SBarry Smith } 2007fabda3fSAlex Fikl shellB->vscale = shellA->vscale; 2017fabda3fSAlex Fikl shellB->vshift = shellA->vshift; 2020c0fd78eSBarry Smith if (shellA->dshift) { 2030c0fd78eSBarry Smith if (!shellB->dshift) { 2040c0fd78eSBarry Smith ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr); 2057fabda3fSAlex Fikl } 2060c0fd78eSBarry Smith ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr); 2077fabda3fSAlex Fikl } else { 2089f137db4SBarry Smith ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr); 2097fabda3fSAlex Fikl } 2100c0fd78eSBarry Smith if (shellA->left) { 2110c0fd78eSBarry Smith if (!shellB->left) { 2120c0fd78eSBarry Smith ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr); 2137fabda3fSAlex Fikl } 2140c0fd78eSBarry Smith ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr); 2157fabda3fSAlex Fikl } else { 2169f137db4SBarry Smith ierr = VecDestroy(&shellB->left);CHKERRQ(ierr); 2177fabda3fSAlex Fikl } 2180c0fd78eSBarry Smith if (shellA->right) { 2190c0fd78eSBarry Smith if (!shellB->right) { 2200c0fd78eSBarry Smith ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr); 2217fabda3fSAlex Fikl } 2220c0fd78eSBarry Smith ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr); 2237fabda3fSAlex Fikl } else { 2249f137db4SBarry Smith ierr = VecDestroy(&shellB->right);CHKERRQ(ierr); 2257fabda3fSAlex Fikl } 22640e381c3SBarry Smith ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr); 22740e381c3SBarry Smith if (shellA->axpy) { 22840e381c3SBarry Smith ierr = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr); 22940e381c3SBarry Smith shellB->axpy = shellA->axpy; 23040e381c3SBarry Smith shellB->axpy_vscale = shellA->axpy_vscale; 23140e381c3SBarry Smith } 2327fabda3fSAlex Fikl PetscFunctionReturn(0); 2337fabda3fSAlex Fikl } 2347fabda3fSAlex Fikl 235cb8c8a77SBarry Smith PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M) 236cb8c8a77SBarry Smith { 237cb8c8a77SBarry Smith PetscErrorCode ierr; 238cb8c8a77SBarry Smith void *ctx; 239cb8c8a77SBarry Smith 240cb8c8a77SBarry Smith PetscFunctionBegin; 241cb8c8a77SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 242cb8c8a77SBarry Smith ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr); 243cb8c8a77SBarry Smith ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 244cb8c8a77SBarry Smith PetscFunctionReturn(0); 245cb8c8a77SBarry Smith } 246cb8c8a77SBarry Smith 247dfbe8321SBarry Smith PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y) 248ef66eb69SBarry Smith { 249ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)A->data; 250dfbe8321SBarry Smith PetscErrorCode ierr; 25125578ef6SJed Brown Vec xx; 252e3f487b0SBarry Smith PetscObjectState instate,outstate; 253ef66eb69SBarry Smith 254ef66eb69SBarry Smith PetscFunctionBegin; 2558fe8eb27SJed Brown ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr); 256e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 257976e8c5aSLisandro Dalcin if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL"); 25828f357e3SAlex Fikl ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr); 259e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 260e3f487b0SBarry Smith if (instate == outstate) { 261e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 262e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 263e3f487b0SBarry Smith } 2648fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 2658fe8eb27SJed Brown ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr); 2669f137db4SBarry Smith 2679f137db4SBarry Smith if (shell->axpy) { 2689f137db4SBarry Smith if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);} 2699f137db4SBarry Smith ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr); 2709f137db4SBarry Smith ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 2719f137db4SBarry Smith } 272ef66eb69SBarry Smith PetscFunctionReturn(0); 273ef66eb69SBarry Smith } 274ef66eb69SBarry Smith 2755edf6516SJed Brown PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z) 2765edf6516SJed Brown { 2775edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 2785edf6516SJed Brown PetscErrorCode ierr; 2795edf6516SJed Brown 2805edf6516SJed Brown PetscFunctionBegin; 2815edf6516SJed Brown if (y == z) { 2825edf6516SJed Brown if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);} 2835edf6516SJed Brown ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr); 284b8430d49SBarry Smith ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr); 2855edf6516SJed Brown } else { 2865edf6516SJed Brown ierr = MatMult(A,x,z);CHKERRQ(ierr); 2875edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 2885edf6516SJed Brown } 2895edf6516SJed Brown PetscFunctionReturn(0); 2905edf6516SJed Brown } 2915edf6516SJed Brown 29218be62a5SSatish Balay PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y) 29318be62a5SSatish Balay { 29418be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 29518be62a5SSatish Balay PetscErrorCode ierr; 29625578ef6SJed Brown Vec xx; 297e3f487b0SBarry Smith PetscObjectState instate,outstate; 29818be62a5SSatish Balay 29918be62a5SSatish Balay PetscFunctionBegin; 3008fe8eb27SJed Brown ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr); 301e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr); 3020c0fd78eSBarry Smith if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL"); 30328f357e3SAlex Fikl ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr); 304e3f487b0SBarry Smith ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr); 305e3f487b0SBarry Smith if (instate == outstate) { 306e3f487b0SBarry Smith /* increase the state of the output vector since the user did not update its state themself as should have been done */ 307e3f487b0SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); 308e3f487b0SBarry Smith } 3098fe8eb27SJed Brown ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr); 3108fe8eb27SJed Brown ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr); 31118be62a5SSatish Balay PetscFunctionReturn(0); 31218be62a5SSatish Balay } 31318be62a5SSatish Balay 3145edf6516SJed Brown PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z) 3155edf6516SJed Brown { 3165edf6516SJed Brown Mat_Shell *shell = (Mat_Shell*)A->data; 3175edf6516SJed Brown PetscErrorCode ierr; 3185edf6516SJed Brown 3195edf6516SJed Brown PetscFunctionBegin; 3205edf6516SJed Brown if (y == z) { 3215edf6516SJed Brown if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);} 3225edf6516SJed Brown ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr); 3235edf6516SJed Brown ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr); 3245edf6516SJed Brown } else { 3255edf6516SJed Brown ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr); 3265edf6516SJed Brown ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr); 3275edf6516SJed Brown } 3285edf6516SJed Brown PetscFunctionReturn(0); 3295edf6516SJed Brown } 3305edf6516SJed Brown 3310c0fd78eSBarry Smith /* 3320c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 3330c0fd78eSBarry Smith */ 33418be62a5SSatish Balay PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v) 33518be62a5SSatish Balay { 33618be62a5SSatish Balay Mat_Shell *shell = (Mat_Shell*)A->data; 33718be62a5SSatish Balay PetscErrorCode ierr; 33818be62a5SSatish Balay 33918be62a5SSatish Balay PetscFunctionBegin; 3400c0fd78eSBarry Smith if (shell->ops->getdiagonal) { 34128f357e3SAlex Fikl ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr); 3420c0fd78eSBarry Smith } else { 3430c0fd78eSBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 3440c0fd78eSBarry Smith } 34518be62a5SSatish Balay ierr = VecScale(v,shell->vscale);CHKERRQ(ierr); 3468fe8eb27SJed Brown if (shell->dshift) { 3470c0fd78eSBarry Smith ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr); 34818be62a5SSatish Balay } 3490c0fd78eSBarry Smith ierr = VecShift(v,shell->vshift);CHKERRQ(ierr); 3508fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);} 3518fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);} 35281c02519SBarry Smith if (shell->axpy) { 35381c02519SBarry Smith if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);} 35481c02519SBarry Smith ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr); 35581c02519SBarry Smith ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr); 35681c02519SBarry Smith } 35718be62a5SSatish Balay PetscFunctionReturn(0); 35818be62a5SSatish Balay } 35918be62a5SSatish Balay 360f4df32b1SMatthew Knepley PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a) 361ef66eb69SBarry Smith { 362ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 3638fe8eb27SJed Brown PetscErrorCode ierr; 364b24ad042SBarry Smith 365ef66eb69SBarry Smith PetscFunctionBegin; 3660c0fd78eSBarry Smith if (shell->left || shell->right) { 3678fe8eb27SJed Brown if (!shell->dshift) { 3680c0fd78eSBarry Smith ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr); 3690c0fd78eSBarry Smith ierr = VecSet(shell->dshift,a);CHKERRQ(ierr); 3700c0fd78eSBarry Smith } else { 3710c0fd78eSBarry Smith if (shell->left) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 3720c0fd78eSBarry Smith if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 3739f137db4SBarry Smith ierr = VecShift(shell->dshift,a);CHKERRQ(ierr); 3740c0fd78eSBarry Smith } 3758fe8eb27SJed Brown if (shell->left) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);} 3768fe8eb27SJed Brown if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);} 3778fe8eb27SJed Brown } else shell->vshift += a; 378ef66eb69SBarry Smith PetscFunctionReturn(0); 379ef66eb69SBarry Smith } 380ef66eb69SBarry Smith 3816464bf51SAlex Fikl PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins) 3826464bf51SAlex Fikl { 3836464bf51SAlex Fikl Mat_Shell *shell = (Mat_Shell*)A->data; 3846464bf51SAlex Fikl PetscErrorCode ierr; 3856464bf51SAlex Fikl 3866464bf51SAlex Fikl PetscFunctionBegin; 3870c0fd78eSBarry Smith if (ins == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported with INSERT_VALUES"); 3880c0fd78eSBarry Smith if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);} 3890c0fd78eSBarry Smith if (shell->left || shell->right) { 39092fabd1fSBarry Smith if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);} 3919f137db4SBarry Smith if (shell->left && shell->right) { 3929f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 3939f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr); 3949f137db4SBarry Smith } else if (shell->left) { 3959f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr); 3969f137db4SBarry Smith } else { 3979f137db4SBarry Smith ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr); 3989f137db4SBarry Smith } 3999f137db4SBarry Smith if (!shell->dshift) { 4009f137db4SBarry Smith ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr); 4019f137db4SBarry Smith ierr = VecCopy(shell->dshift,shell->right_work);CHKERRQ(ierr); 4029f137db4SBarry Smith } else { 4039f137db4SBarry Smith ierr = VecAXPY(shell->dshift,1.0,shell->right_work);CHKERRQ(ierr); 4049f137db4SBarry Smith } 4050c0fd78eSBarry Smith } else { 4060c0fd78eSBarry Smith ierr = VecAXPY(shell->dshift,1.0,D);CHKERRQ(ierr); 4076464bf51SAlex Fikl } 4086464bf51SAlex Fikl PetscFunctionReturn(0); 4096464bf51SAlex Fikl } 4106464bf51SAlex Fikl 411f4df32b1SMatthew Knepley PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a) 412ef66eb69SBarry Smith { 413ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 4148fe8eb27SJed Brown PetscErrorCode ierr; 415b24ad042SBarry Smith 416ef66eb69SBarry Smith PetscFunctionBegin; 417f4df32b1SMatthew Knepley shell->vscale *= a; 4180c0fd78eSBarry Smith shell->vshift *= a; 4192205254eSKarl Rupp if (shell->dshift) { 4202205254eSKarl Rupp ierr = VecScale(shell->dshift,a);CHKERRQ(ierr); 4210c0fd78eSBarry Smith } 42281c02519SBarry Smith shell->axpy_vscale *= a; 4238fe8eb27SJed Brown PetscFunctionReturn(0); 42418be62a5SSatish Balay } 4258fe8eb27SJed Brown 4268fe8eb27SJed Brown static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right) 4278fe8eb27SJed Brown { 4288fe8eb27SJed Brown Mat_Shell *shell = (Mat_Shell*)Y->data; 4298fe8eb27SJed Brown PetscErrorCode ierr; 4308fe8eb27SJed Brown 4318fe8eb27SJed Brown PetscFunctionBegin; 43281c02519SBarry Smith if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation"); 4338fe8eb27SJed Brown if (left) { 4340c0fd78eSBarry Smith if (!shell->left) { 4350c0fd78eSBarry Smith ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr); 4368fe8eb27SJed Brown ierr = VecCopy(left,shell->left);CHKERRQ(ierr); 4370c0fd78eSBarry Smith } else { 4380c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr); 43918be62a5SSatish Balay } 440ef66eb69SBarry Smith } 4418fe8eb27SJed Brown if (right) { 4420c0fd78eSBarry Smith if (!shell->right) { 4430c0fd78eSBarry Smith ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr); 4448fe8eb27SJed Brown ierr = VecCopy(right,shell->right);CHKERRQ(ierr); 4450c0fd78eSBarry Smith } else { 4460c0fd78eSBarry Smith ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr); 4478fe8eb27SJed Brown } 4488fe8eb27SJed Brown } 449ef66eb69SBarry Smith PetscFunctionReturn(0); 450ef66eb69SBarry Smith } 451ef66eb69SBarry Smith 452dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t) 453ef66eb69SBarry Smith { 454ef66eb69SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 4550c0fd78eSBarry Smith PetscErrorCode ierr; 456ef66eb69SBarry Smith 457ef66eb69SBarry Smith PetscFunctionBegin; 4588fe8eb27SJed Brown if (t == MAT_FINAL_ASSEMBLY) { 459ef66eb69SBarry Smith shell->vshift = 0.0; 460ef66eb69SBarry Smith shell->vscale = 1.0; 4610c0fd78eSBarry Smith ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr); 4620c0fd78eSBarry Smith ierr = VecDestroy(&shell->left);CHKERRQ(ierr); 4630c0fd78eSBarry Smith ierr = VecDestroy(&shell->right);CHKERRQ(ierr); 46481c02519SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 465ef66eb69SBarry Smith } 466ef66eb69SBarry Smith PetscFunctionReturn(0); 467ef66eb69SBarry Smith } 468ef66eb69SBarry Smith 469cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*); 470b951964fSBarry Smith 4713b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool *missing,PetscInt *d) 4723b49f96aSBarry Smith { 4733b49f96aSBarry Smith PetscFunctionBegin; 4743b49f96aSBarry Smith *missing = PETSC_FALSE; 4753b49f96aSBarry Smith PetscFunctionReturn(0); 4763b49f96aSBarry Smith } 4773b49f96aSBarry Smith 4789f137db4SBarry Smith PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str) 4799f137db4SBarry Smith { 4809f137db4SBarry Smith Mat_Shell *shell = (Mat_Shell*)Y->data; 4819f137db4SBarry Smith PetscErrorCode ierr; 4829f137db4SBarry Smith 4839f137db4SBarry Smith PetscFunctionBegin; 4849f137db4SBarry Smith ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 4859f137db4SBarry Smith ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr); 4869f137db4SBarry Smith shell->axpy = X; 4879f137db4SBarry Smith shell->axpy_vscale = a; 4889f137db4SBarry Smith PetscFunctionReturn(0); 4899f137db4SBarry Smith } 4909f137db4SBarry Smith 49109dc0095SBarry Smith static struct _MatOps MatOps_Values = {0, 49220563c6bSBarry Smith 0, 49320563c6bSBarry Smith 0, 4949f137db4SBarry Smith 0, 4950c0fd78eSBarry Smith /* 4*/ MatMultAdd_Shell, 4969f137db4SBarry Smith 0, 4970c0fd78eSBarry Smith MatMultTransposeAdd_Shell, 498b951964fSBarry Smith 0, 499b951964fSBarry Smith 0, 500b951964fSBarry Smith 0, 50197304618SKris Buschelman /*10*/ 0, 502b951964fSBarry Smith 0, 503b951964fSBarry Smith 0, 504b951964fSBarry Smith 0, 505b951964fSBarry Smith 0, 50697304618SKris Buschelman /*15*/ 0, 507b951964fSBarry Smith 0, 5080c0fd78eSBarry Smith MatGetDiagonal_Shell, 5098fe8eb27SJed Brown MatDiagonalScale_Shell, 510b951964fSBarry Smith 0, 51197304618SKris Buschelman /*20*/ 0, 512ef66eb69SBarry Smith MatAssemblyEnd_Shell, 513b951964fSBarry Smith 0, 514b951964fSBarry Smith 0, 515d519adbfSMatthew Knepley /*24*/ 0, 516b951964fSBarry Smith 0, 517b951964fSBarry Smith 0, 518b951964fSBarry Smith 0, 519b951964fSBarry Smith 0, 520d519adbfSMatthew Knepley /*29*/ 0, 521b951964fSBarry Smith 0, 522273d9f13SBarry Smith 0, 523b951964fSBarry Smith 0, 524b951964fSBarry Smith 0, 525cb8c8a77SBarry Smith /*34*/ MatDuplicate_Shell, 526b951964fSBarry Smith 0, 527b951964fSBarry Smith 0, 52809dc0095SBarry Smith 0, 52909dc0095SBarry Smith 0, 5309f137db4SBarry Smith /*39*/ MatAXPY_Shell, 53109dc0095SBarry Smith 0, 53209dc0095SBarry Smith 0, 53309dc0095SBarry Smith 0, 534cb8c8a77SBarry Smith MatCopy_Shell, 535d519adbfSMatthew Knepley /*44*/ 0, 536ef66eb69SBarry Smith MatScale_Shell, 537ef66eb69SBarry Smith MatShift_Shell, 5386464bf51SAlex Fikl MatDiagonalSet_Shell, 53909dc0095SBarry Smith 0, 540f73d5cc4SBarry Smith /*49*/ 0, 54109dc0095SBarry Smith 0, 54209dc0095SBarry Smith 0, 54309dc0095SBarry Smith 0, 54409dc0095SBarry Smith 0, 545d519adbfSMatthew Knepley /*54*/ 0, 54609dc0095SBarry Smith 0, 54709dc0095SBarry Smith 0, 54809dc0095SBarry Smith 0, 54909dc0095SBarry Smith 0, 550d519adbfSMatthew Knepley /*59*/ 0, 551b9b97703SBarry Smith MatDestroy_Shell, 55209dc0095SBarry Smith 0, 553357abbc8SBarry Smith 0, 554273d9f13SBarry Smith 0, 555d519adbfSMatthew Knepley /*64*/ 0, 556273d9f13SBarry Smith 0, 557273d9f13SBarry Smith 0, 558273d9f13SBarry Smith 0, 559273d9f13SBarry Smith 0, 560d519adbfSMatthew Knepley /*69*/ 0, 561273d9f13SBarry Smith 0, 562c87e5d42SMatthew Knepley MatConvert_Shell, 563273d9f13SBarry Smith 0, 56497304618SKris Buschelman 0, 565d519adbfSMatthew Knepley /*74*/ 0, 56697304618SKris Buschelman 0, 56797304618SKris Buschelman 0, 56897304618SKris Buschelman 0, 56997304618SKris Buschelman 0, 570d519adbfSMatthew Knepley /*79*/ 0, 57197304618SKris Buschelman 0, 57297304618SKris Buschelman 0, 57397304618SKris Buschelman 0, 574865e5f61SKris Buschelman 0, 575d519adbfSMatthew Knepley /*84*/ 0, 576865e5f61SKris Buschelman 0, 577865e5f61SKris Buschelman 0, 578865e5f61SKris Buschelman 0, 579865e5f61SKris Buschelman 0, 580d519adbfSMatthew Knepley /*89*/ 0, 581865e5f61SKris Buschelman 0, 582865e5f61SKris Buschelman 0, 583865e5f61SKris Buschelman 0, 584865e5f61SKris Buschelman 0, 585d519adbfSMatthew Knepley /*94*/ 0, 586865e5f61SKris Buschelman 0, 587865e5f61SKris Buschelman 0, 5883964eb88SJed Brown 0, 5893964eb88SJed Brown 0, 5903964eb88SJed Brown /*99*/ 0, 5913964eb88SJed Brown 0, 5923964eb88SJed Brown 0, 5933964eb88SJed Brown 0, 5943964eb88SJed Brown 0, 5953964eb88SJed Brown /*104*/ 0, 5963964eb88SJed Brown 0, 5973964eb88SJed Brown 0, 5983964eb88SJed Brown 0, 5993964eb88SJed Brown 0, 6003964eb88SJed Brown /*109*/ 0, 6013964eb88SJed Brown 0, 6023964eb88SJed Brown 0, 6033964eb88SJed Brown 0, 6043b49f96aSBarry Smith MatMissingDiagonal_Shell, 6053964eb88SJed Brown /*114*/ 0, 6063964eb88SJed Brown 0, 6073964eb88SJed Brown 0, 6083964eb88SJed Brown 0, 6093964eb88SJed Brown 0, 6103964eb88SJed Brown /*119*/ 0, 6113964eb88SJed Brown 0, 6123964eb88SJed Brown 0, 6133964eb88SJed Brown 0, 6143964eb88SJed Brown 0, 6153964eb88SJed Brown /*124*/ 0, 6163964eb88SJed Brown 0, 6173964eb88SJed Brown 0, 6183964eb88SJed Brown 0, 6193964eb88SJed Brown 0, 6203964eb88SJed Brown /*129*/ 0, 6213964eb88SJed Brown 0, 6223964eb88SJed Brown 0, 6233964eb88SJed Brown 0, 6243964eb88SJed Brown 0, 6253964eb88SJed Brown /*134*/ 0, 6263964eb88SJed Brown 0, 6273964eb88SJed Brown 0, 6283964eb88SJed Brown 0, 6293964eb88SJed Brown 0, 6303964eb88SJed Brown /*139*/ 0, 631f9426fe0SMark Adams 0, 6323964eb88SJed Brown 0 6333964eb88SJed Brown }; 634273d9f13SBarry Smith 6350bad9183SKris Buschelman /*MC 636fafad747SKris Buschelman MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free. 6370bad9183SKris Buschelman 6380bad9183SKris Buschelman Level: advanced 6390bad9183SKris Buschelman 6400c0fd78eSBarry Smith .seealso: MatCreateShell() 6410bad9183SKris Buschelman M*/ 6420bad9183SKris Buschelman 6438cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A) 644273d9f13SBarry Smith { 645273d9f13SBarry Smith Mat_Shell *b; 646dfbe8321SBarry Smith PetscErrorCode ierr; 647273d9f13SBarry Smith 648273d9f13SBarry Smith PetscFunctionBegin; 649273d9f13SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 650273d9f13SBarry Smith 651b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 652273d9f13SBarry Smith A->data = (void*)b; 653273d9f13SBarry Smith 65426283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 65526283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 656273d9f13SBarry Smith 657273d9f13SBarry Smith b->ctx = 0; 658ef66eb69SBarry Smith b->vshift = 0.0; 659ef66eb69SBarry Smith b->vscale = 1.0; 6600c0fd78eSBarry Smith b->managescalingshifts = PETSC_TRUE; 661273d9f13SBarry Smith A->assembled = PETSC_TRUE; 662210f0121SBarry Smith A->preallocated = PETSC_FALSE; 6632205254eSKarl Rupp 66417667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr); 665273d9f13SBarry Smith PetscFunctionReturn(0); 666273d9f13SBarry Smith } 667e51e0e81SBarry Smith 6684b828684SBarry Smith /*@C 669052efed2SBarry Smith MatCreateShell - Creates a new matrix class for use with a user-defined 670ff756334SLois Curfman McInnes private data storage format. 671e51e0e81SBarry Smith 672c7fcc2eaSBarry Smith Collective on MPI_Comm 673c7fcc2eaSBarry Smith 674e51e0e81SBarry Smith Input Parameters: 675c7fcc2eaSBarry Smith + comm - MPI communicator 676c7fcc2eaSBarry Smith . m - number of local rows (must be given) 677c7fcc2eaSBarry Smith . n - number of local columns (must be given) 678c7fcc2eaSBarry Smith . M - number of global rows (may be PETSC_DETERMINE) 679c7fcc2eaSBarry Smith . N - number of global columns (may be PETSC_DETERMINE) 680c7fcc2eaSBarry Smith - ctx - pointer to data needed by the shell matrix routines 681e51e0e81SBarry Smith 682ff756334SLois Curfman McInnes Output Parameter: 68344cd7ae7SLois Curfman McInnes . A - the matrix 684e51e0e81SBarry Smith 685ff2fd236SBarry Smith Level: advanced 686ff2fd236SBarry Smith 687f39d1f56SLois Curfman McInnes Usage: 6887b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 689f39d1f56SLois Curfman McInnes $ MatCreateShell(comm,m,n,M,N,ctx,&mat); 690c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 691f39d1f56SLois Curfman McInnes $ [ Use matrix for operations that have been set ] 692f39d1f56SLois Curfman McInnes $ MatDestroy(mat); 693f39d1f56SLois Curfman McInnes 694ff756334SLois Curfman McInnes Notes: 695ff756334SLois Curfman McInnes The shell matrix type is intended to provide a simple class to use 696ff756334SLois Curfman McInnes with KSP (such as, for use with matrix-free methods). You should not 697ff756334SLois Curfman McInnes use the shell type if you plan to define a complete matrix class. 698e51e0e81SBarry Smith 699*95452b02SPatrick Sanan Fortran Notes: 700*95452b02SPatrick Sanan To use this from Fortran with a ctx you must write an interface definition for this 701daf670e6SBarry Smith function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing 702daf670e6SBarry Smith in as the ctx argument. 703f38a8d56SBarry Smith 704f39d1f56SLois Curfman McInnes PETSc requires that matrices and vectors being used for certain 705f39d1f56SLois Curfman McInnes operations are partitioned accordingly. For example, when 706645985a0SLois Curfman McInnes creating a shell matrix, A, that supports parallel matrix-vector 707645985a0SLois Curfman McInnes products using MatMult(A,x,y) the user should set the number 708645985a0SLois Curfman McInnes of local matrix rows to be the number of local elements of the 709645985a0SLois Curfman McInnes corresponding result vector, y. Note that this is information is 710645985a0SLois Curfman McInnes required for use of the matrix interface routines, even though 711645985a0SLois Curfman McInnes the shell matrix may not actually be physically partitioned. 712645985a0SLois Curfman McInnes For example, 713f39d1f56SLois Curfman McInnes 714f39d1f56SLois Curfman McInnes $ 715f39d1f56SLois Curfman McInnes $ Vec x, y 7167b2a1423SBarry Smith $ extern int mult(Mat,Vec,Vec); 717645985a0SLois Curfman McInnes $ Mat A 718f39d1f56SLois Curfman McInnes $ 719c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,M,&y); 720c94f878dSBarry Smith $ VecCreateMPI(comm,PETSC_DECIDE,N,&x); 721f39d1f56SLois Curfman McInnes $ VecGetLocalSize(y,&m); 722c7fcc2eaSBarry Smith $ VecGetLocalSize(x,&n); 723c7fcc2eaSBarry Smith $ MatCreateShell(comm,m,n,M,N,ctx,&A); 724c134de8dSSatish Balay $ MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult); 725645985a0SLois Curfman McInnes $ MatMult(A,x,y); 726645985a0SLois Curfman McInnes $ MatDestroy(A); 727f39d1f56SLois Curfman McInnes $ VecDestroy(y); VecDestroy(x); 728645985a0SLois Curfman McInnes $ 729e51e0e81SBarry Smith 7309b53d762SBarry Smith 7319b53d762SBarry Smith MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), and MatScale() internally so these 7329b53d762SBarry Smith operations cannot be overwritten unless MatShellSetManageScalingShifts() is called. 7339b53d762SBarry Smith 7349b53d762SBarry Smith 7350c0fd78eSBarry Smith For rectangular matrices do all the scalings and shifts make sense? 7360c0fd78eSBarry Smith 737*95452b02SPatrick Sanan Developers Notes: 738*95452b02SPatrick Sanan Regarding shifting and scaling. The general form is 7390c0fd78eSBarry Smith 7400c0fd78eSBarry Smith diag(left)(vscale*A + diag(dshift) + vshift I)diag(right) 7410c0fd78eSBarry Smith 7420c0fd78eSBarry Smith The order you apply the operations is important. For example if you have a dshift then 7430c0fd78eSBarry Smith apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift 7440c0fd78eSBarry Smith you get s*vscale*A + diag(shift) 7450c0fd78eSBarry Smith 7460c0fd78eSBarry Smith A is the user provided function. 7470c0fd78eSBarry Smith 7480b627109SLois Curfman McInnes .keywords: matrix, shell, create 7490b627109SLois Curfman McInnes 7500c0fd78eSBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts() 751e51e0e81SBarry Smith @*/ 7527087cfbeSBarry Smith PetscErrorCode MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A) 753e51e0e81SBarry Smith { 754dfbe8321SBarry Smith PetscErrorCode ierr; 755ed3cc1f0SBarry Smith 7563a40ed3dSBarry Smith PetscFunctionBegin; 757f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 758f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 759273d9f13SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 760273d9f13SBarry Smith ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr); 7610fcd0fccSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 762273d9f13SBarry Smith PetscFunctionReturn(0); 763c7fcc2eaSBarry Smith } 764c7fcc2eaSBarry Smith 765c6866cfdSSatish Balay /*@ 766273d9f13SBarry Smith MatShellSetContext - sets the context for a shell matrix 767c7fcc2eaSBarry Smith 7683f9fe445SBarry Smith Logically Collective on Mat 769c7fcc2eaSBarry Smith 770273d9f13SBarry Smith Input Parameters: 771273d9f13SBarry Smith + mat - the shell matrix 772273d9f13SBarry Smith - ctx - the context 773273d9f13SBarry Smith 774273d9f13SBarry Smith Level: advanced 775273d9f13SBarry Smith 776*95452b02SPatrick Sanan Fortran Notes: 777*95452b02SPatrick Sanan To use this from Fortran you must write a Fortran interface definition for this 778daf670e6SBarry Smith function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 779273d9f13SBarry Smith 780273d9f13SBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation() 7810bc0a719SBarry Smith @*/ 7827087cfbeSBarry Smith PetscErrorCode MatShellSetContext(Mat mat,void *ctx) 783273d9f13SBarry Smith { 784273d9f13SBarry Smith Mat_Shell *shell = (Mat_Shell*)mat->data; 785dfbe8321SBarry Smith PetscErrorCode ierr; 786ace3abfcSBarry Smith PetscBool flg; 787273d9f13SBarry Smith 788273d9f13SBarry Smith PetscFunctionBegin; 7890700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 790251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 791273d9f13SBarry Smith if (flg) { 792273d9f13SBarry Smith shell->ctx = ctx; 793ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix"); 7943a40ed3dSBarry Smith PetscFunctionReturn(0); 795e51e0e81SBarry Smith } 796e51e0e81SBarry Smith 7970c0fd78eSBarry Smith /*@ 7980c0fd78eSBarry Smith MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately 7990c0fd78eSBarry Smith after MatCreateShell() 8000c0fd78eSBarry Smith 8010c0fd78eSBarry Smith Logically Collective on Mat 8020c0fd78eSBarry Smith 8030c0fd78eSBarry Smith Input Parameter: 8040c0fd78eSBarry Smith . mat - the shell matrix 8050c0fd78eSBarry Smith 8060c0fd78eSBarry Smith Level: advanced 8070c0fd78eSBarry Smith 8080c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation() 8090c0fd78eSBarry Smith @*/ 8100c0fd78eSBarry Smith PetscErrorCode MatShellSetManageScalingShifts(Mat A) 8110c0fd78eSBarry Smith { 8120c0fd78eSBarry Smith PetscErrorCode ierr; 813976e8c5aSLisandro Dalcin Mat_Shell *shell; 8140c0fd78eSBarry Smith PetscBool flg; 8150c0fd78eSBarry Smith 8160c0fd78eSBarry Smith PetscFunctionBegin; 8170c0fd78eSBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 8180c0fd78eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr); 8190c0fd78eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 820976e8c5aSLisandro Dalcin shell = (Mat_Shell*)A->data; 8210c0fd78eSBarry Smith shell->managescalingshifts = PETSC_FALSE; 822976e8c5aSLisandro Dalcin A->ops->diagonalset = NULL; 823976e8c5aSLisandro Dalcin A->ops->diagonalscale = NULL; 824976e8c5aSLisandro Dalcin A->ops->scale = NULL; 825976e8c5aSLisandro Dalcin A->ops->shift = NULL; 826976e8c5aSLisandro Dalcin A->ops->axpy = NULL; 8270c0fd78eSBarry Smith PetscFunctionReturn(0); 8280c0fd78eSBarry Smith } 8290c0fd78eSBarry Smith 830c16cb8f2SBarry Smith /*@C 831f3b1f45cSBarry Smith MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function. 832f3b1f45cSBarry Smith 833f3b1f45cSBarry Smith Logically Collective on Mat 834f3b1f45cSBarry Smith 835f3b1f45cSBarry Smith Input Parameters: 836f3b1f45cSBarry Smith + mat - the shell matrix 837f3b1f45cSBarry Smith . f - the function 838f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 839f3b1f45cSBarry Smith - ctx - an optional context for the function 840f3b1f45cSBarry Smith 841f3b1f45cSBarry Smith Output Parameter: 842f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 843f3b1f45cSBarry Smith 844f3b1f45cSBarry Smith Options Database: 845f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 846f3b1f45cSBarry Smith 847f3b1f45cSBarry Smith Level: advanced 848f3b1f45cSBarry Smith 849*95452b02SPatrick Sanan Fortran Notes: 850*95452b02SPatrick Sanan Not supported from Fortran 851f3b1f45cSBarry Smith 852f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose() 853f3b1f45cSBarry Smith @*/ 854f3b1f45cSBarry Smith PetscErrorCode MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 855f3b1f45cSBarry Smith { 856f3b1f45cSBarry Smith PetscErrorCode ierr; 857f3b1f45cSBarry Smith PetscInt m,n; 858f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 859f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 86074e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 861f3b1f45cSBarry Smith 862f3b1f45cSBarry Smith PetscFunctionBegin; 863f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 864f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr); 865f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 866f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr); 867f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 868f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 869f3b1f45cSBarry Smith 870f3b1f45cSBarry Smith ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr); 871f3b1f45cSBarry Smith ierr = MatComputeExplicitOperator(mat,&Dmat);CHKERRQ(ierr); 872f3b1f45cSBarry Smith 873f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 874f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 875f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 876f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 877f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 878f3b1f45cSBarry Smith flag = PETSC_FALSE; 879f3b1f45cSBarry Smith if (v) { 880fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm));CHKERRQ(ierr); 881f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 882f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 883f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr); 884f3b1f45cSBarry Smith } 885f3b1f45cSBarry Smith } else if (v) { 886fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr); 887f3b1f45cSBarry Smith } 888f3b1f45cSBarry Smith if (flg) *flg = flag; 889f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 890f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 891f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 892f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 893f3b1f45cSBarry Smith PetscFunctionReturn(0); 894f3b1f45cSBarry Smith } 895f3b1f45cSBarry Smith 896f3b1f45cSBarry Smith /*@C 897f3b1f45cSBarry Smith MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function. 898f3b1f45cSBarry Smith 899f3b1f45cSBarry Smith Logically Collective on Mat 900f3b1f45cSBarry Smith 901f3b1f45cSBarry Smith Input Parameters: 902f3b1f45cSBarry Smith + mat - the shell matrix 903f3b1f45cSBarry Smith . f - the function 904f3b1f45cSBarry Smith . base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated 905f3b1f45cSBarry Smith - ctx - an optional context for the function 906f3b1f45cSBarry Smith 907f3b1f45cSBarry Smith Output Parameter: 908f3b1f45cSBarry Smith . flg - PETSC_TRUE if the multiply is likely correct 909f3b1f45cSBarry Smith 910f3b1f45cSBarry Smith Options Database: 911f3b1f45cSBarry Smith . -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference 912f3b1f45cSBarry Smith 913f3b1f45cSBarry Smith Level: advanced 914f3b1f45cSBarry Smith 915*95452b02SPatrick Sanan Fortran Notes: 916*95452b02SPatrick Sanan Not supported from Fortran 917f3b1f45cSBarry Smith 918f3b1f45cSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult() 919f3b1f45cSBarry Smith @*/ 920f3b1f45cSBarry Smith PetscErrorCode MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg) 921f3b1f45cSBarry Smith { 922f3b1f45cSBarry Smith PetscErrorCode ierr; 923f3b1f45cSBarry Smith Vec x,y,z; 924f3b1f45cSBarry Smith PetscInt m,n,M,N; 925f3b1f45cSBarry Smith Mat mf,Dmf,Dmat,Ddiff; 926f3b1f45cSBarry Smith PetscReal Diffnorm,Dmfnorm; 92774e5cdcaSBarry Smith PetscBool v = PETSC_FALSE, flag = PETSC_TRUE; 928f3b1f45cSBarry Smith 929f3b1f45cSBarry Smith PetscFunctionBegin; 930f3b1f45cSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 931f3b1f45cSBarry Smith ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr); 932f3b1f45cSBarry Smith ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr); 933f3b1f45cSBarry Smith ierr = VecDuplicate(y,&z);CHKERRQ(ierr); 934f3b1f45cSBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 935f3b1f45cSBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 936f3b1f45cSBarry Smith ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr); 937f3b1f45cSBarry Smith ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr); 938f3b1f45cSBarry Smith ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr); 939f3b1f45cSBarry Smith ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr); 940f3b1f45cSBarry Smith ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr); 941f3b1f45cSBarry Smith ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr); 942f3b1f45cSBarry Smith 943f3b1f45cSBarry Smith ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr); 944f3b1f45cSBarry Smith ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 945f3b1f45cSBarry Smith ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr); 946f3b1f45cSBarry Smith ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr); 947f3b1f45cSBarry Smith if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) { 948f3b1f45cSBarry Smith flag = PETSC_FALSE; 949f3b1f45cSBarry Smith if (v) { 950fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm));CHKERRQ(ierr); 951f3b1f45cSBarry Smith ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 952f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 953f3b1f45cSBarry Smith ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr); 954f3b1f45cSBarry Smith } 955f3b1f45cSBarry Smith } else if (v) { 956fc7aafd1SBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr); 957f3b1f45cSBarry Smith } 958f3b1f45cSBarry Smith if (flg) *flg = flag; 959f3b1f45cSBarry Smith ierr = MatDestroy(&mf);CHKERRQ(ierr); 960f3b1f45cSBarry Smith ierr = MatDestroy(&Dmat);CHKERRQ(ierr); 961f3b1f45cSBarry Smith ierr = MatDestroy(&Ddiff);CHKERRQ(ierr); 962f3b1f45cSBarry Smith ierr = MatDestroy(&Dmf);CHKERRQ(ierr); 963f3b1f45cSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 964f3b1f45cSBarry Smith ierr = VecDestroy(&y);CHKERRQ(ierr); 965f3b1f45cSBarry Smith ierr = VecDestroy(&z);CHKERRQ(ierr); 966f3b1f45cSBarry Smith PetscFunctionReturn(0); 967f3b1f45cSBarry Smith } 968f3b1f45cSBarry Smith 969f3b1f45cSBarry Smith /*@C 9700c0fd78eSBarry Smith MatShellSetOperation - Allows user to set a matrix operation for a shell matrix. 971e51e0e81SBarry Smith 9723f9fe445SBarry Smith Logically Collective on Mat 973fee21e36SBarry Smith 974c7fcc2eaSBarry Smith Input Parameters: 975c7fcc2eaSBarry Smith + mat - the shell matrix 976c7fcc2eaSBarry Smith . op - the name of the operation 977c7fcc2eaSBarry Smith - f - the function that provides the operation. 978c7fcc2eaSBarry Smith 97915091d37SBarry Smith Level: advanced 98015091d37SBarry Smith 981fae171e0SBarry Smith Usage: 982e93bc3c1Svictor $ extern PetscErrorCode usermult(Mat,Vec,Vec); 983f39d1f56SLois Curfman McInnes $ ierr = MatCreateShell(comm,m,n,M,N,ctx,&A); 984c134de8dSSatish Balay $ ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult); 9850b627109SLois Curfman McInnes 986a62d957aSLois Curfman McInnes Notes: 987e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 9881c1c02c0SLois Curfman McInnes operations, which all have the form MATOP_<OPERATION>, where 989a62d957aSLois Curfman McInnes <OPERATION> is the name (in all capital letters) of the 9901c1c02c0SLois Curfman McInnes user interface routine (e.g., MatMult() -> MATOP_MULT). 991a62d957aSLois Curfman McInnes 992a4d85fd7SBarry Smith All user-provided functions (except for MATOP_DESTROY) should have the same calling 993deebb3c3SLois Curfman McInnes sequence as the usual matrix interface routines, since they 994deebb3c3SLois Curfman McInnes are intended to be accessed via the usual matrix interface 995deebb3c3SLois Curfman McInnes routines, e.g., 996a62d957aSLois Curfman McInnes $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 997a62d957aSLois Curfman McInnes 9984aa34b0aSBarry Smith In particular each function MUST return an error code of 0 on success and 9994aa34b0aSBarry Smith nonzero on failure. 10004aa34b0aSBarry Smith 1001a62d957aSLois Curfman McInnes Within each user-defined routine, the user should call 1002a62d957aSLois Curfman McInnes MatShellGetContext() to obtain the user-defined context that was 1003a62d957aSLois Curfman McInnes set by MatCreateShell(). 1004a62d957aSLois Curfman McInnes 1005*95452b02SPatrick Sanan Fortran Notes: 1006*95452b02SPatrick Sanan For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not 1007501d9185SBarry Smith generate a matrix. See src/mat/examples/tests/ex120f.F 1008501d9185SBarry Smith 10090c0fd78eSBarry Smith Use MatSetOperation() to set an operation for any matrix type 10100c0fd78eSBarry Smith 1011a62d957aSLois Curfman McInnes .keywords: matrix, shell, set, operation 1012a62d957aSLois Curfman McInnes 10130c0fd78eSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts() 1014e51e0e81SBarry Smith @*/ 10157087cfbeSBarry Smith PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void)) 1016e51e0e81SBarry Smith { 1017ace3abfcSBarry Smith PetscBool flg; 1018976e8c5aSLisandro Dalcin Mat_Shell *shell; 1019976e8c5aSLisandro Dalcin PetscErrorCode ierr; 1020273d9f13SBarry Smith 10213a40ed3dSBarry Smith PetscFunctionBegin; 10220700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 10230c0fd78eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 10240c0fd78eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 1025976e8c5aSLisandro Dalcin shell = (Mat_Shell*)mat->data; 1026976e8c5aSLisandro Dalcin 10275edf6516SJed Brown switch (op) { 10285edf6516SJed Brown case MATOP_DESTROY: 102928f357e3SAlex Fikl shell->ops->destroy = (PetscErrorCode (*)(Mat))f; 10305edf6516SJed Brown break; 1031976e8c5aSLisandro Dalcin case MATOP_VIEW: 1032976e8c5aSLisandro Dalcin if (!mat->ops->viewnative) { 1033976e8c5aSLisandro Dalcin mat->ops->viewnative = mat->ops->view; 1034976e8c5aSLisandro Dalcin } 1035976e8c5aSLisandro Dalcin mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f; 1036976e8c5aSLisandro Dalcin break; 1037976e8c5aSLisandro Dalcin case MATOP_COPY: 1038976e8c5aSLisandro Dalcin shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f; 1039976e8c5aSLisandro Dalcin break; 10406464bf51SAlex Fikl case MATOP_DIAGONAL_SET: 10410c0fd78eSBarry Smith case MATOP_DIAGONAL_SCALE: 10420c0fd78eSBarry Smith case MATOP_SHIFT: 10430c0fd78eSBarry Smith case MATOP_SCALE: 10449f137db4SBarry Smith case MATOP_AXPY: 10450c0fd78eSBarry Smith if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()"); 10460c0fd78eSBarry Smith (((void(**)(void))mat->ops)[op]) = f; 10476464bf51SAlex Fikl break; 10480c0fd78eSBarry Smith case MATOP_GET_DIAGONAL: 1049976e8c5aSLisandro Dalcin if (shell->managescalingshifts) { 1050976e8c5aSLisandro Dalcin shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 1051976e8c5aSLisandro Dalcin mat->ops->getdiagonal = MatGetDiagonal_Shell; 1052976e8c5aSLisandro Dalcin } else { 1053976e8c5aSLisandro Dalcin shell->ops->getdiagonal = NULL; 1054976e8c5aSLisandro Dalcin mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f; 105540a5a6c4SBarry Smith } 10565edf6516SJed Brown break; 10575edf6516SJed Brown case MATOP_MULT: 10589f137db4SBarry Smith if (shell->managescalingshifts) { 10599f137db4SBarry Smith shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 10609f137db4SBarry Smith mat->ops->mult = MatMult_Shell; 1061976e8c5aSLisandro Dalcin } else { 1062976e8c5aSLisandro Dalcin shell->ops->mult = NULL; 1063976e8c5aSLisandro Dalcin mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1064976e8c5aSLisandro Dalcin } 10655edf6516SJed Brown break; 10665edf6516SJed Brown case MATOP_MULT_TRANSPOSE: 10679f137db4SBarry Smith if (shell->managescalingshifts) { 10689f137db4SBarry Smith shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 10695259c5a4SBarry Smith mat->ops->multtranspose = MatMultTranspose_Shell; 1070976e8c5aSLisandro Dalcin } else { 1071976e8c5aSLisandro Dalcin shell->ops->multtranspose = NULL; 1072976e8c5aSLisandro Dalcin mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f; 1073976e8c5aSLisandro Dalcin } 10745edf6516SJed Brown break; 10755edf6516SJed Brown default: 10765edf6516SJed Brown (((void(**)(void))mat->ops)[op]) = f; 10775ab264f3SAlp Dener break; 1078a62d957aSLois Curfman McInnes } 10793a40ed3dSBarry Smith PetscFunctionReturn(0); 1080e51e0e81SBarry Smith } 1081f0479e8cSBarry Smith 1082d4bb536fSBarry Smith /*@C 1083d4bb536fSBarry Smith MatShellGetOperation - Gets a matrix function for a shell matrix. 1084d4bb536fSBarry Smith 1085c7fcc2eaSBarry Smith Not Collective 1086c7fcc2eaSBarry Smith 1087d4bb536fSBarry Smith Input Parameters: 1088c7fcc2eaSBarry Smith + mat - the shell matrix 1089c7fcc2eaSBarry Smith - op - the name of the operation 1090d4bb536fSBarry Smith 1091d4bb536fSBarry Smith Output Parameter: 1092d4bb536fSBarry Smith . f - the function that provides the operation. 1093d4bb536fSBarry Smith 109415091d37SBarry Smith Level: advanced 109515091d37SBarry Smith 1096d4bb536fSBarry Smith Notes: 1097e090d566SSatish Balay See the file include/petscmat.h for a complete list of matrix 1098d4bb536fSBarry Smith operations, which all have the form MATOP_<OPERATION>, where 1099d4bb536fSBarry Smith <OPERATION> is the name (in all capital letters) of the 1100d4bb536fSBarry Smith user interface routine (e.g., MatMult() -> MATOP_MULT). 1101d4bb536fSBarry Smith 1102d4bb536fSBarry Smith All user-provided functions have the same calling 1103d4bb536fSBarry Smith sequence as the usual matrix interface routines, since they 1104d4bb536fSBarry Smith are intended to be accessed via the usual matrix interface 1105d4bb536fSBarry Smith routines, e.g., 1106d4bb536fSBarry Smith $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec) 1107d4bb536fSBarry Smith 1108d4bb536fSBarry Smith Within each user-defined routine, the user should call 1109d4bb536fSBarry Smith MatShellGetContext() to obtain the user-defined context that was 1110d4bb536fSBarry Smith set by MatCreateShell(). 1111d4bb536fSBarry Smith 1112d4bb536fSBarry Smith .keywords: matrix, shell, set, operation 1113d4bb536fSBarry Smith 1114ab50ec6bSBarry Smith .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext() 1115d4bb536fSBarry Smith @*/ 11167087cfbeSBarry Smith PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void)) 1117d4bb536fSBarry Smith { 1118ace3abfcSBarry Smith PetscBool flg; 1119976e8c5aSLisandro Dalcin Mat_Shell *shell; 1120976e8c5aSLisandro Dalcin PetscErrorCode ierr; 1121273d9f13SBarry Smith 11223a40ed3dSBarry Smith PetscFunctionBegin; 11230700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1124976e8c5aSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr); 1125976e8c5aSLisandro Dalcin if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices"); 1126976e8c5aSLisandro Dalcin shell = (Mat_Shell*)mat->data; 1127976e8c5aSLisandro Dalcin 112828f357e3SAlex Fikl switch (op) { 112928f357e3SAlex Fikl case MATOP_DESTROY: 113028f357e3SAlex Fikl *f = (void (*)(void))shell->ops->destroy; 113128f357e3SAlex Fikl break; 113228f357e3SAlex Fikl case MATOP_VIEW: 113328f357e3SAlex Fikl *f = (void (*)(void))mat->ops->view; 113428f357e3SAlex Fikl break; 1135976e8c5aSLisandro Dalcin case MATOP_COPY: 1136976e8c5aSLisandro Dalcin *f = (void (*)(void))shell->ops->copy; 1137976e8c5aSLisandro Dalcin break; 1138976e8c5aSLisandro Dalcin case MATOP_DIAGONAL_SET: 1139976e8c5aSLisandro Dalcin case MATOP_DIAGONAL_SCALE: 1140976e8c5aSLisandro Dalcin case MATOP_SHIFT: 1141976e8c5aSLisandro Dalcin case MATOP_SCALE: 1142976e8c5aSLisandro Dalcin case MATOP_AXPY: 1143976e8c5aSLisandro Dalcin *f = (((void (**)(void))mat->ops)[op]); 1144976e8c5aSLisandro Dalcin break; 1145976e8c5aSLisandro Dalcin case MATOP_GET_DIAGONAL: 1146976e8c5aSLisandro Dalcin if (shell->ops->getdiagonal) 1147976e8c5aSLisandro Dalcin *f = (void (*)(void))shell->ops->getdiagonal; 1148976e8c5aSLisandro Dalcin else 1149976e8c5aSLisandro Dalcin *f = (((void (**)(void))mat->ops)[op]); 1150976e8c5aSLisandro Dalcin break; 1151976e8c5aSLisandro Dalcin case MATOP_MULT: 1152976e8c5aSLisandro Dalcin if (shell->ops->mult) 1153976e8c5aSLisandro Dalcin *f = (void (*)(void))shell->ops->mult; 1154976e8c5aSLisandro Dalcin else 1155976e8c5aSLisandro Dalcin *f = (((void (**)(void))mat->ops)[op]); 1156976e8c5aSLisandro Dalcin break; 1157976e8c5aSLisandro Dalcin case MATOP_MULT_TRANSPOSE: 1158976e8c5aSLisandro Dalcin if (shell->ops->multtranspose) 1159976e8c5aSLisandro Dalcin *f = (void (*)(void))shell->ops->multtranspose; 1160976e8c5aSLisandro Dalcin else 1161976e8c5aSLisandro Dalcin *f = (((void (**)(void))mat->ops)[op]); 1162976e8c5aSLisandro Dalcin break; 116328f357e3SAlex Fikl default: 1164c134de8dSSatish Balay *f = (((void (**)(void))mat->ops)[op]); 1165d4bb536fSBarry Smith } 11663a40ed3dSBarry Smith PetscFunctionReturn(0); 1167d4bb536fSBarry Smith } 1168