1*ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 2aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h> /*I "submatfree.h" I*/ 3a7e14dcfSSatish Balay 4a7e14dcfSSatish Balay #undef __FUNCT__ 5a7e14dcfSSatish Balay #define __FUNCT__ "MatCreateSubMatrixFree" 6a7e14dcfSSatish Balay /*@C 7a7e14dcfSSatish Balay MatCreateSubMatrixFree - Creates a reduced matrix by masking a 8a7e14dcfSSatish Balay full matrix. 9a7e14dcfSSatish Balay 10a7e14dcfSSatish Balay Collective on matrix 11a7e14dcfSSatish Balay 12a7e14dcfSSatish Balay Input Parameters: 13a7e14dcfSSatish Balay + mat - matrix of arbitrary type 14a7e14dcfSSatish Balay . Rows - the rows that will be in the submatrix 15a7e14dcfSSatish Balay - Cols - the columns that will be in the submatrix 16a7e14dcfSSatish Balay 17a7e14dcfSSatish Balay Output Parameters: 18a7e14dcfSSatish Balay . J - New matrix 19a7e14dcfSSatish Balay 20a7e14dcfSSatish Balay Notes: 21a7e14dcfSSatish Balay The user provides the input data and is responsible for destroying 22a7e14dcfSSatish Balay this data after matrix J has been destroyed. 23a7e14dcfSSatish Balay 24a7e14dcfSSatish Balay Level: developer 25a7e14dcfSSatish Balay 26a7e14dcfSSatish Balay .seealso: MatCreate() 27a7e14dcfSSatish Balay @*/ 28a7e14dcfSSatish Balay PetscErrorCode MatCreateSubMatrixFree(Mat mat,IS Rows, IS Cols, Mat *J) 29a7e14dcfSSatish Balay { 30a7e14dcfSSatish Balay MPI_Comm comm=((PetscObject)mat)->comm; 31a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 32a7e14dcfSSatish Balay PetscErrorCode ierr; 33a7e14dcfSSatish Balay PetscMPIInt size; 34a7e14dcfSSatish Balay PetscInt mloc,nloc,m,n; 35a7e14dcfSSatish Balay 3687f595a5SBarry Smith PetscFunctionBegin; 373c9e27cfSGeoffrey Irving ierr = PetscNew(&ctx);CHKERRQ(ierr); 38a7e14dcfSSatish Balay ctx->A=mat; 39a7e14dcfSSatish Balay ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 40a7e14dcfSSatish Balay ierr = MatGetLocalSize(mat,&mloc,&nloc);CHKERRQ(ierr); 41a7e14dcfSSatish Balay ierr = MPI_Comm_size(comm,&size); 42a7e14dcfSSatish Balay if (size == 1) { 43a7e14dcfSSatish Balay ierr = VecCreateSeq(comm,n,&ctx->VC);CHKERRQ(ierr); 44a7e14dcfSSatish Balay } else { 45a7e14dcfSSatish Balay ierr = VecCreateMPI(comm,nloc,n,&ctx->VC);CHKERRQ(ierr); 46a7e14dcfSSatish Balay } 47a7e14dcfSSatish Balay ctx->VR=ctx->VC; 48a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 49a7e14dcfSSatish Balay 50a7e14dcfSSatish Balay 51a7e14dcfSSatish Balay ctx->Rows = Rows; 52a7e14dcfSSatish Balay ctx->Cols = Cols; 53a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)Rows);CHKERRQ(ierr); 54a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)Cols);CHKERRQ(ierr); 55a7e14dcfSSatish Balay ierr = MatCreateShell(comm,mloc,nloc,m,n,ctx,J);CHKERRQ(ierr); 56a7e14dcfSSatish Balay 57a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))MatMult_SMF);CHKERRQ(ierr); 58a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))MatDestroy_SMF);CHKERRQ(ierr); 59a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))MatView_SMF);CHKERRQ(ierr); 60a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_SMF);CHKERRQ(ierr); 61a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DIAGONAL_SET,(void(*)(void))MatDiagonalSet_SMF);CHKERRQ(ierr); 62a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_SHIFT,(void(*)(void))MatShift_SMF);CHKERRQ(ierr); 63a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_EQUAL,(void(*)(void))MatEqual_SMF);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_SCALE,(void(*)(void))MatScale_SMF);CHKERRQ(ierr); 65a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_TRANSPOSE,(void(*)(void))MatTranspose_SMF);CHKERRQ(ierr); 66a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_SMF);CHKERRQ(ierr); 67a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_GET_SUBMATRICES,(void(*)(void))MatGetSubMatrices_SMF);CHKERRQ(ierr); 68a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_NORM,(void(*)(void))MatNorm_SMF);CHKERRQ(ierr); 69a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_SMF);CHKERRQ(ierr); 70a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_GET_SUBMATRIX,(void(*)(void))MatGetSubMatrix_SMF);CHKERRQ(ierr); 71a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_GET_ROW_MAX,(void(*)(void))MatDuplicate_SMF);CHKERRQ(ierr); 72a7e14dcfSSatish Balay 73a7e14dcfSSatish Balay ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)(*J));CHKERRQ(ierr); 74a7e14dcfSSatish Balay PetscFunctionReturn(0); 75a7e14dcfSSatish Balay } 76a7e14dcfSSatish Balay 77a7e14dcfSSatish Balay #undef __FUNCT__ 78a7e14dcfSSatish Balay #define __FUNCT__ "MatSMFResetRowColumn" 79a7e14dcfSSatish Balay PetscErrorCode MatSMFResetRowColumn(Mat mat,IS Rows,IS Cols){ 80a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 81a7e14dcfSSatish Balay PetscErrorCode ierr; 82a7e14dcfSSatish Balay 83a7e14dcfSSatish Balay PetscFunctionBegin; 84a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 85a7e14dcfSSatish Balay ierr = ISDestroy(&ctx->Rows);CHKERRQ(ierr); 86a7e14dcfSSatish Balay ierr = ISDestroy(&ctx->Cols);CHKERRQ(ierr); 87a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)Rows);CHKERRQ(ierr); 88a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)Cols);CHKERRQ(ierr); 89a7e14dcfSSatish Balay ctx->Cols=Cols; 90a7e14dcfSSatish Balay ctx->Rows=Rows; 91a7e14dcfSSatish Balay PetscFunctionReturn(0); 92a7e14dcfSSatish Balay } 93a7e14dcfSSatish Balay 94a7e14dcfSSatish Balay #undef __FUNCT__ 95a7e14dcfSSatish Balay #define __FUNCT__ "MatMult_SMF" 96a7e14dcfSSatish Balay PetscErrorCode MatMult_SMF(Mat mat,Vec a,Vec y) 97a7e14dcfSSatish Balay { 98a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 99a7e14dcfSSatish Balay PetscErrorCode ierr; 100a7e14dcfSSatish Balay 101a7e14dcfSSatish Balay PetscFunctionBegin; 102a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 103a7e14dcfSSatish Balay ierr = VecCopy(a,ctx->VR);CHKERRQ(ierr); 1044473680cSBarry Smith ierr = VecISSet(ctx->VR,ctx->Cols,0.0);CHKERRQ(ierr); 105a7e14dcfSSatish Balay ierr = MatMult(ctx->A,ctx->VR,y);CHKERRQ(ierr); 1064473680cSBarry Smith ierr = VecISSet(y,ctx->Rows,0.0);CHKERRQ(ierr); 107a7e14dcfSSatish Balay PetscFunctionReturn(0); 108a7e14dcfSSatish Balay } 109a7e14dcfSSatish Balay 110a7e14dcfSSatish Balay #undef __FUNCT__ 111a7e14dcfSSatish Balay #define __FUNCT__ "MatMultTranspose_SMF" 112a7e14dcfSSatish Balay PetscErrorCode MatMultTranspose_SMF(Mat mat,Vec a,Vec y) 113a7e14dcfSSatish Balay { 114a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 115a7e14dcfSSatish Balay PetscErrorCode ierr; 116a7e14dcfSSatish Balay 117a7e14dcfSSatish Balay PetscFunctionBegin; 118a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 119a7e14dcfSSatish Balay ierr = VecCopy(a,ctx->VC);CHKERRQ(ierr); 1204473680cSBarry Smith ierr = VecISSet(ctx->VC,ctx->Rows,0.0);CHKERRQ(ierr); 121a7e14dcfSSatish Balay ierr = MatMultTranspose(ctx->A,ctx->VC,y);CHKERRQ(ierr); 1224473680cSBarry Smith ierr = VecISSet(y,ctx->Cols,0.0);CHKERRQ(ierr); 123a7e14dcfSSatish Balay PetscFunctionReturn(0); 124a7e14dcfSSatish Balay } 125a7e14dcfSSatish Balay 126a7e14dcfSSatish Balay #undef __FUNCT__ 127a7e14dcfSSatish Balay #define __FUNCT__ "MatDiagonalSet_SMF" 128a7e14dcfSSatish Balay PetscErrorCode MatDiagonalSet_SMF(Mat M, Vec D,InsertMode is) 129a7e14dcfSSatish Balay { 130a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 131a7e14dcfSSatish Balay PetscErrorCode ierr; 132a7e14dcfSSatish Balay 133a7e14dcfSSatish Balay PetscFunctionBegin; 134a7e14dcfSSatish Balay ierr = MatShellGetContext(M,(void **)&ctx);CHKERRQ(ierr); 135a7e14dcfSSatish Balay ierr = MatDiagonalSet(ctx->A,D,is); 136a7e14dcfSSatish Balay PetscFunctionReturn(0); 137a7e14dcfSSatish Balay } 138a7e14dcfSSatish Balay 139a7e14dcfSSatish Balay #undef __FUNCT__ 140a7e14dcfSSatish Balay #define __FUNCT__ "MatDestroy_SMF" 141a7e14dcfSSatish Balay PetscErrorCode MatDestroy_SMF(Mat mat) 142a7e14dcfSSatish Balay { 143a7e14dcfSSatish Balay PetscErrorCode ierr; 144a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 145a7e14dcfSSatish Balay 146a7e14dcfSSatish Balay PetscFunctionBegin; 147a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 148a7e14dcfSSatish Balay ierr = MatDestroy(&ctx->A);CHKERRQ(ierr); 149a7e14dcfSSatish Balay ierr = ISDestroy(&ctx->Rows);CHKERRQ(ierr); 150a7e14dcfSSatish Balay ierr = ISDestroy(&ctx->Cols);CHKERRQ(ierr); 151a7e14dcfSSatish Balay ierr = VecDestroy(&ctx->VC);CHKERRQ(ierr); 152a7e14dcfSSatish Balay ierr = PetscFree(ctx);CHKERRQ(ierr); 153a7e14dcfSSatish Balay PetscFunctionReturn(0); 154a7e14dcfSSatish Balay } 155a7e14dcfSSatish Balay 156a7e14dcfSSatish Balay 157a7e14dcfSSatish Balay 158a7e14dcfSSatish Balay #undef __FUNCT__ 159a7e14dcfSSatish Balay #define __FUNCT__ "MatView_SMF" 160a7e14dcfSSatish Balay PetscErrorCode MatView_SMF(Mat mat,PetscViewer viewer) 161a7e14dcfSSatish Balay { 162a7e14dcfSSatish Balay PetscErrorCode ierr; 163a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 164a7e14dcfSSatish Balay 165a7e14dcfSSatish Balay PetscFunctionBegin; 166a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 167a7e14dcfSSatish Balay ierr = MatView(ctx->A,viewer);CHKERRQ(ierr); 168a7e14dcfSSatish Balay PetscFunctionReturn(0); 169a7e14dcfSSatish Balay } 170a7e14dcfSSatish Balay 171a7e14dcfSSatish Balay #undef __FUNCT__ 172a7e14dcfSSatish Balay #define __FUNCT__ "MatShift_SMF" 173a7e14dcfSSatish Balay PetscErrorCode MatShift_SMF(Mat Y, PetscReal a) 174a7e14dcfSSatish Balay { 175a7e14dcfSSatish Balay PetscErrorCode ierr; 176a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 177a7e14dcfSSatish Balay 178a7e14dcfSSatish Balay PetscFunctionBegin; 179a7e14dcfSSatish Balay ierr = MatShellGetContext(Y,(void **)&ctx);CHKERRQ(ierr); 180a7e14dcfSSatish Balay ierr = MatShift(ctx->A,a);CHKERRQ(ierr); 181a7e14dcfSSatish Balay PetscFunctionReturn(0); 182a7e14dcfSSatish Balay } 183a7e14dcfSSatish Balay 184a7e14dcfSSatish Balay #undef __FUNCT__ 185a7e14dcfSSatish Balay #define __FUNCT__ "MatDuplicate_SMF" 186a7e14dcfSSatish Balay PetscErrorCode MatDuplicate_SMF(Mat mat,MatDuplicateOption op,Mat *M) 187a7e14dcfSSatish Balay { 188a7e14dcfSSatish Balay PetscErrorCode ierr; 189a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 190a7e14dcfSSatish Balay 191a7e14dcfSSatish Balay PetscFunctionBegin; 192a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 193a7e14dcfSSatish Balay ierr = MatCreateSubMatrixFree(ctx->A,ctx->Rows,ctx->Cols,M);CHKERRQ(ierr); 194a7e14dcfSSatish Balay PetscFunctionReturn(0); 195a7e14dcfSSatish Balay } 196a7e14dcfSSatish Balay 197a7e14dcfSSatish Balay #undef __FUNCT__ 198a7e14dcfSSatish Balay #define __FUNCT__ "MatEqual_SMF" 199a7e14dcfSSatish Balay PetscErrorCode MatEqual_SMF(Mat A,Mat B,PetscBool *flg) 200a7e14dcfSSatish Balay { 201a7e14dcfSSatish Balay PetscErrorCode ierr; 202a7e14dcfSSatish Balay MatSubMatFreeCtx ctx1,ctx2; 203a7e14dcfSSatish Balay PetscBool flg1,flg2,flg3; 204a7e14dcfSSatish Balay 205a7e14dcfSSatish Balay PetscFunctionBegin; 206a7e14dcfSSatish Balay ierr = MatShellGetContext(A,(void **)&ctx1);CHKERRQ(ierr); 207a7e14dcfSSatish Balay ierr = MatShellGetContext(B,(void **)&ctx2);CHKERRQ(ierr); 208a7e14dcfSSatish Balay ierr = ISEqual(ctx1->Rows,ctx2->Rows,&flg2);CHKERRQ(ierr); 209a7e14dcfSSatish Balay ierr = ISEqual(ctx1->Cols,ctx2->Cols,&flg3);CHKERRQ(ierr); 210a7e14dcfSSatish Balay if (flg2==PETSC_FALSE || flg3==PETSC_FALSE){ 211a7e14dcfSSatish Balay *flg=PETSC_FALSE; 212a7e14dcfSSatish Balay } else { 213a7e14dcfSSatish Balay ierr = MatEqual(ctx1->A,ctx2->A,&flg1);CHKERRQ(ierr); 214a7e14dcfSSatish Balay if (flg1==PETSC_FALSE){ *flg=PETSC_FALSE;} 215a7e14dcfSSatish Balay else { *flg=PETSC_TRUE;} 216a7e14dcfSSatish Balay } 217a7e14dcfSSatish Balay PetscFunctionReturn(0); 218a7e14dcfSSatish Balay } 219a7e14dcfSSatish Balay 220a7e14dcfSSatish Balay #undef __FUNCT__ 221a7e14dcfSSatish Balay #define __FUNCT__ "MatScale_SMF" 222a7e14dcfSSatish Balay PetscErrorCode MatScale_SMF(Mat mat, PetscReal a) 223a7e14dcfSSatish Balay { 224a7e14dcfSSatish Balay PetscErrorCode ierr; 225a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 226a7e14dcfSSatish Balay 227a7e14dcfSSatish Balay PetscFunctionBegin; 228a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 229a7e14dcfSSatish Balay ierr = MatScale(ctx->A,a);CHKERRQ(ierr); 230a7e14dcfSSatish Balay PetscFunctionReturn(0); 231a7e14dcfSSatish Balay } 232a7e14dcfSSatish Balay 233a7e14dcfSSatish Balay #undef __FUNCT__ 234a7e14dcfSSatish Balay #define __FUNCT__ "MatTranspose_SMF" 235a7e14dcfSSatish Balay PetscErrorCode MatTranspose_SMF(Mat mat,Mat *B) 236a7e14dcfSSatish Balay { 237a7e14dcfSSatish Balay PetscFunctionBegin; 238a7e14dcfSSatish Balay PetscFunctionReturn(1); 239a7e14dcfSSatish Balay } 240a7e14dcfSSatish Balay 241a7e14dcfSSatish Balay #undef __FUNCT__ 242a7e14dcfSSatish Balay #define __FUNCT__ "MatGetDiagonal_SMF" 243a7e14dcfSSatish Balay PetscErrorCode MatGetDiagonal_SMF(Mat mat,Vec v) 244a7e14dcfSSatish Balay { 245a7e14dcfSSatish Balay PetscErrorCode ierr; 246a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 247a7e14dcfSSatish Balay 248a7e14dcfSSatish Balay PetscFunctionBegin; 249a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 250a7e14dcfSSatish Balay ierr = MatGetDiagonal(ctx->A,v);CHKERRQ(ierr); 251a7e14dcfSSatish Balay PetscFunctionReturn(0); 252a7e14dcfSSatish Balay } 253a7e14dcfSSatish Balay 254a7e14dcfSSatish Balay #undef __FUNCT__ 255a7e14dcfSSatish Balay #define __FUNCT__ "MatDiagonalSet_SMF" 256a7e14dcfSSatish Balay PetscErrorCode MatGetRowMax_SMF(Mat M, Vec D) 257a7e14dcfSSatish Balay { 258a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 259a7e14dcfSSatish Balay PetscErrorCode ierr; 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay PetscFunctionBegin; 262a7e14dcfSSatish Balay ierr = MatShellGetContext(M,(void **)&ctx);CHKERRQ(ierr); 2636c23d075SBarry Smith ierr = MatGetRowMax(ctx->A,D,NULL); 264a7e14dcfSSatish Balay PetscFunctionReturn(0); 265a7e14dcfSSatish Balay } 266a7e14dcfSSatish Balay 267a7e14dcfSSatish Balay #undef __FUNCT__ 268a7e14dcfSSatish Balay #define __FUNCT__ "MatGetSubMatrices_SMF" 269a7e14dcfSSatish Balay PetscErrorCode MatGetSubMatrices_SMF(Mat A,PetscInt n, IS *irow,IS *icol,MatReuse scall,Mat **B) 270a7e14dcfSSatish Balay { 271a7e14dcfSSatish Balay PetscErrorCode ierr; 272a7e14dcfSSatish Balay PetscInt i; 273a7e14dcfSSatish Balay 274a7e14dcfSSatish Balay PetscFunctionBegin; 275a7e14dcfSSatish Balay if (scall == MAT_INITIAL_MATRIX) { 2760e660641SBarry Smith ierr = PetscMalloc1(n+1,B );CHKERRQ(ierr); 277a7e14dcfSSatish Balay } 278a7e14dcfSSatish Balay 279a7e14dcfSSatish Balay for ( i=0; i<n; i++ ) { 280a7e14dcfSSatish Balay ierr = MatGetSubMatrix_SMF(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 281a7e14dcfSSatish Balay } 282a7e14dcfSSatish Balay PetscFunctionReturn(0); 283a7e14dcfSSatish Balay } 284a7e14dcfSSatish Balay 285a7e14dcfSSatish Balay #undef __FUNCT__ 286a7e14dcfSSatish Balay #define __FUNCT__ "MatGetSubMatrix_SMF" 287a7e14dcfSSatish Balay PetscErrorCode MatGetSubMatrix_SMF(Mat mat,IS isrow,IS iscol,MatReuse cll, 288a7e14dcfSSatish Balay Mat *newmat) 289a7e14dcfSSatish Balay { 290a7e14dcfSSatish Balay PetscErrorCode ierr; 291a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 292a7e14dcfSSatish Balay 293a7e14dcfSSatish Balay PetscFunctionBegin; 294a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 295a7e14dcfSSatish Balay if (newmat){ 296a7e14dcfSSatish Balay ierr = MatDestroy(&*newmat);CHKERRQ(ierr); 297a7e14dcfSSatish Balay } 298a7e14dcfSSatish Balay ierr = MatCreateSubMatrixFree(ctx->A,isrow,iscol, newmat);CHKERRQ(ierr); 299a7e14dcfSSatish Balay PetscFunctionReturn(0); 300a7e14dcfSSatish Balay } 301a7e14dcfSSatish Balay 302a7e14dcfSSatish Balay #undef __FUNCT__ 303a7e14dcfSSatish Balay #define __FUNCT__ "MatGetRow_SMF" 304a7e14dcfSSatish Balay PetscErrorCode MatGetRow_SMF(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt **cols,const PetscReal **vals) 305a7e14dcfSSatish Balay { 306a7e14dcfSSatish Balay PetscErrorCode ierr; 307a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 308a7e14dcfSSatish Balay 309a7e14dcfSSatish Balay PetscFunctionBegin; 310a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 311a7e14dcfSSatish Balay ierr = MatGetRow(ctx->A,row,ncols,cols,vals);CHKERRQ(ierr); 312a7e14dcfSSatish Balay PetscFunctionReturn(0); 313a7e14dcfSSatish Balay } 314a7e14dcfSSatish Balay 315a7e14dcfSSatish Balay #undef __FUNCT__ 316a7e14dcfSSatish Balay #define __FUNCT__ "MatRestoreRow_SMF" 317a7e14dcfSSatish Balay PetscErrorCode MatRestoreRow_SMF(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt **cols,const PetscReal **vals) 318a7e14dcfSSatish Balay { 319a7e14dcfSSatish Balay PetscErrorCode ierr; 320a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 321a7e14dcfSSatish Balay 322a7e14dcfSSatish Balay PetscFunctionBegin; 323a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 324a7e14dcfSSatish Balay ierr = MatRestoreRow(ctx->A,row,ncols,cols,vals);CHKERRQ(ierr); 325a7e14dcfSSatish Balay PetscFunctionReturn(0); 326a7e14dcfSSatish Balay } 327a7e14dcfSSatish Balay 328a7e14dcfSSatish Balay #undef __FUNCT__ 329a7e14dcfSSatish Balay #define __FUNCT__ "MatGetColumnVector_SMF" 330a7e14dcfSSatish Balay PetscErrorCode MatGetColumnVector_SMF(Mat mat,Vec Y, PetscInt col) 331a7e14dcfSSatish Balay { 332a7e14dcfSSatish Balay PetscErrorCode ierr; 333a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 334a7e14dcfSSatish Balay 335a7e14dcfSSatish Balay PetscFunctionBegin; 336a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 337a7e14dcfSSatish Balay ierr = MatGetColumnVector(ctx->A,Y,col);CHKERRQ(ierr); 338a7e14dcfSSatish Balay PetscFunctionReturn(0); 339a7e14dcfSSatish Balay } 340a7e14dcfSSatish Balay 341a7e14dcfSSatish Balay #undef __FUNCT__ 342a7e14dcfSSatish Balay #define __FUNCT__ "MatConvert_SMF" 343a7e14dcfSSatish Balay PetscErrorCode MatConvert_SMF(Mat mat,MatType newtype,Mat *NewMat) 344a7e14dcfSSatish Balay { 345a7e14dcfSSatish Balay PetscErrorCode ierr; 346a7e14dcfSSatish Balay PetscMPIInt size; 347a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 348a7e14dcfSSatish Balay 349a7e14dcfSSatish Balay PetscFunctionBegin; 350a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 35187f595a5SBarry Smith ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr); 352a7e14dcfSSatish Balay PetscFunctionReturn(1); 353a7e14dcfSSatish Balay } 354a7e14dcfSSatish Balay 355a7e14dcfSSatish Balay #undef __FUNCT__ 356a7e14dcfSSatish Balay #define __FUNCT__ "MatNorm_SMF" 357a7e14dcfSSatish Balay PetscErrorCode MatNorm_SMF(Mat mat,NormType type,PetscReal *norm) 358a7e14dcfSSatish Balay { 359a7e14dcfSSatish Balay PetscErrorCode ierr; 360a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 361a7e14dcfSSatish Balay 362a7e14dcfSSatish Balay PetscFunctionBegin; 363a7e14dcfSSatish Balay ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 364a7e14dcfSSatish Balay if (type == NORM_FROBENIUS) { 365a7e14dcfSSatish Balay *norm = 1.0; 366a7e14dcfSSatish Balay } else if (type == NORM_1 || type == NORM_INFINITY) { 367a7e14dcfSSatish Balay *norm = 1.0; 36887f595a5SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 369a7e14dcfSSatish Balay PetscFunctionReturn(0); 370a7e14dcfSSatish Balay } 371a7e14dcfSSatish Balay 372