1ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 2954e39ddSJose E. Roman #include <../src/tao/matrix/submatfree.h> 3a7e14dcfSSatish Balay 4cc4c1da9SBarry Smith /*@ 5a7e14dcfSSatish Balay MatCreateSubMatrixFree - Creates a reduced matrix by masking a 6a7e14dcfSSatish Balay full matrix. 7a7e14dcfSSatish Balay 8c3339decSBarry Smith Collective 9a7e14dcfSSatish Balay 10a7e14dcfSSatish Balay Input Parameters: 11a7e14dcfSSatish Balay + mat - matrix of arbitrary type 12a7e14dcfSSatish Balay . Rows - the rows that will be in the submatrix 13a7e14dcfSSatish Balay - Cols - the columns that will be in the submatrix 14a7e14dcfSSatish Balay 152fe279fdSBarry Smith Output Parameter: 16a7e14dcfSSatish Balay . J - New matrix 17a7e14dcfSSatish Balay 18a7e14dcfSSatish Balay Level: developer 19a7e14dcfSSatish Balay 20c3339decSBarry Smith Note: 21c3339decSBarry Smith The caller is responsible for destroying the input objects after matrix J has been destroyed. 22c3339decSBarry Smith 23a3b724e8SBarry Smith Developer Note: 24a3b724e8SBarry Smith This should be moved/supported in `Mat` 25a3b724e8SBarry Smith 26db781477SPatrick Sanan .seealso: `MatCreate()` 27a7e14dcfSSatish Balay @*/ 28d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrixFree(Mat mat, IS Rows, IS Cols, Mat *J) 29d71ae5a4SJacob Faibussowitsch { 30367daffbSBarry Smith MPI_Comm comm = PetscObjectComm((PetscObject)mat); 31a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 32a7e14dcfSSatish Balay PetscInt mloc, nloc, m, n; 33a7e14dcfSSatish Balay 3487f595a5SBarry Smith PetscFunctionBegin; 359566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 36a7e14dcfSSatish Balay ctx->A = mat; 379566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat, &m, &n)); 389566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat, &mloc, &nloc)); 399566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat, NULL, &ctx->VC)); 40a7e14dcfSSatish Balay ctx->VR = ctx->VC; 419566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 42a7e14dcfSSatish Balay 43a7e14dcfSSatish Balay ctx->Rows = Rows; 44a7e14dcfSSatish Balay ctx->Cols = Cols; 459566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Rows)); 469566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Cols)); 479566063dSJacob Faibussowitsch PetscCall(MatCreateShell(comm, mloc, nloc, m, n, ctx, J)); 489566063dSJacob Faibussowitsch PetscCall(MatShellSetManageScalingShifts(*J)); 49*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_MULT, (PetscErrorCodeFn *)MatMult_SMF)); 50*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (PetscErrorCodeFn *)MatDestroy_SMF)); 51*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_VIEW, (PetscErrorCodeFn *)MatView_SMF)); 52*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_SMF)); 53*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_DIAGONAL_SET, (PetscErrorCodeFn *)MatDiagonalSet_SMF)); 54*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_SHIFT, (PetscErrorCodeFn *)MatShift_SMF)); 55*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_EQUAL, (PetscErrorCodeFn *)MatEqual_SMF)); 56*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_SCALE, (PetscErrorCodeFn *)MatScale_SMF)); 57*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_TRANSPOSE, (PetscErrorCodeFn *)MatTranspose_SMF)); 58*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_SMF)); 59*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRICES, (PetscErrorCodeFn *)MatCreateSubMatrices_SMF)); 60*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_NORM, (PetscErrorCodeFn *)MatNorm_SMF)); 61*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_DUPLICATE, (PetscErrorCodeFn *)MatDuplicate_SMF)); 62*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRIX, (PetscErrorCodeFn *)MatCreateSubMatrix_SMF)); 63*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_GET_ROW_MAX, (PetscErrorCodeFn *)MatDuplicate_SMF)); 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65a7e14dcfSSatish Balay } 66a7e14dcfSSatish Balay 67d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSMFResetRowColumn(Mat mat, IS Rows, IS Cols) 68d71ae5a4SJacob Faibussowitsch { 69a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 70a7e14dcfSSatish Balay 71a7e14dcfSSatish Balay PetscFunctionBegin; 729566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx->Rows)); 749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx->Cols)); 759566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Rows)); 769566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Cols)); 77a7e14dcfSSatish Balay ctx->Cols = Cols; 78a7e14dcfSSatish Balay ctx->Rows = Rows; 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80a7e14dcfSSatish Balay } 81a7e14dcfSSatish Balay 82d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_SMF(Mat mat, Vec a, Vec y) 83d71ae5a4SJacob Faibussowitsch { 84a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 85a7e14dcfSSatish Balay 86a7e14dcfSSatish Balay PetscFunctionBegin; 879566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 889566063dSJacob Faibussowitsch PetscCall(VecCopy(a, ctx->VR)); 899566063dSJacob Faibussowitsch PetscCall(VecISSet(ctx->VR, ctx->Cols, 0.0)); 909566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->A, ctx->VR, y)); 919566063dSJacob Faibussowitsch PetscCall(VecISSet(y, ctx->Rows, 0.0)); 923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 93a7e14dcfSSatish Balay } 94a7e14dcfSSatish Balay 95d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTranspose_SMF(Mat mat, Vec a, Vec y) 96d71ae5a4SJacob Faibussowitsch { 97a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 98a7e14dcfSSatish Balay 99a7e14dcfSSatish Balay PetscFunctionBegin; 1009566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 1019566063dSJacob Faibussowitsch PetscCall(VecCopy(a, ctx->VC)); 1029566063dSJacob Faibussowitsch PetscCall(VecISSet(ctx->VC, ctx->Rows, 0.0)); 1039566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ctx->A, ctx->VC, y)); 1049566063dSJacob Faibussowitsch PetscCall(VecISSet(y, ctx->Cols, 0.0)); 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106a7e14dcfSSatish Balay } 107a7e14dcfSSatish Balay 108d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet_SMF(Mat M, Vec D, InsertMode is) 109d71ae5a4SJacob Faibussowitsch { 110a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 111a7e14dcfSSatish Balay 112a7e14dcfSSatish Balay PetscFunctionBegin; 1139566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(M, &ctx)); 1149566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(ctx->A, D, is)); 1153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116a7e14dcfSSatish Balay } 117a7e14dcfSSatish Balay 118d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_SMF(Mat mat) 119d71ae5a4SJacob Faibussowitsch { 120a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 121a7e14dcfSSatish Balay 122a7e14dcfSSatish Balay PetscFunctionBegin; 1239566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 1249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->A)); 1259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx->Rows)); 1269566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx->Cols)); 1279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->VC)); 1289566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 130a7e14dcfSSatish Balay } 131a7e14dcfSSatish Balay 132d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SMF(Mat mat, PetscViewer viewer) 133d71ae5a4SJacob Faibussowitsch { 134a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay PetscFunctionBegin; 1379566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 1389566063dSJacob Faibussowitsch PetscCall(MatView(ctx->A, viewer)); 1393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 140a7e14dcfSSatish Balay } 141a7e14dcfSSatish Balay 142d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift_SMF(Mat Y, PetscReal a) 143d71ae5a4SJacob Faibussowitsch { 144a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 145a7e14dcfSSatish Balay 146a7e14dcfSSatish Balay PetscFunctionBegin; 1479566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(Y, &ctx)); 1489566063dSJacob Faibussowitsch PetscCall(MatShift(ctx->A, a)); 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150a7e14dcfSSatish Balay } 151a7e14dcfSSatish Balay 152d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_SMF(Mat mat, MatDuplicateOption op, Mat *M) 153d71ae5a4SJacob Faibussowitsch { 154a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 155a7e14dcfSSatish Balay 156a7e14dcfSSatish Balay PetscFunctionBegin; 1579566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 1589566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrixFree(ctx->A, ctx->Rows, ctx->Cols, M)); 1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 160a7e14dcfSSatish Balay } 161a7e14dcfSSatish Balay 162d71ae5a4SJacob Faibussowitsch PetscErrorCode MatEqual_SMF(Mat A, Mat B, PetscBool *flg) 163d71ae5a4SJacob Faibussowitsch { 164a7e14dcfSSatish Balay MatSubMatFreeCtx ctx1, ctx2; 165a7e14dcfSSatish Balay PetscBool flg1, flg2, flg3; 166a7e14dcfSSatish Balay 167a7e14dcfSSatish Balay PetscFunctionBegin; 1689566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &ctx1)); 1699566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(B, &ctx2)); 1709566063dSJacob Faibussowitsch PetscCall(ISEqual(ctx1->Rows, ctx2->Rows, &flg2)); 1719566063dSJacob Faibussowitsch PetscCall(ISEqual(ctx1->Cols, ctx2->Cols, &flg3)); 172a7e14dcfSSatish Balay if (flg2 == PETSC_FALSE || flg3 == PETSC_FALSE) { 173a7e14dcfSSatish Balay *flg = PETSC_FALSE; 174a7e14dcfSSatish Balay } else { 1759566063dSJacob Faibussowitsch PetscCall(MatEqual(ctx1->A, ctx2->A, &flg1)); 1769371c9d4SSatish Balay if (flg1 == PETSC_FALSE) { 1779371c9d4SSatish Balay *flg = PETSC_FALSE; 1789371c9d4SSatish Balay } else { 1799371c9d4SSatish Balay *flg = PETSC_TRUE; 1809371c9d4SSatish Balay } 181a7e14dcfSSatish Balay } 1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 183a7e14dcfSSatish Balay } 184a7e14dcfSSatish Balay 185d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_SMF(Mat mat, PetscReal a) 186d71ae5a4SJacob Faibussowitsch { 187a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 188a7e14dcfSSatish Balay 189a7e14dcfSSatish Balay PetscFunctionBegin; 1909566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 1919566063dSJacob Faibussowitsch PetscCall(MatScale(ctx->A, a)); 1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 193a7e14dcfSSatish Balay } 194a7e14dcfSSatish Balay 195d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTranspose_SMF(Mat mat, Mat *B) 196d71ae5a4SJacob Faibussowitsch { 197a7e14dcfSSatish Balay PetscFunctionBegin; 19811cc89d2SBarry Smith SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "No support for transpose for MatCreateSubMatrixFree() matrix"); 1993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 200a7e14dcfSSatish Balay } 201a7e14dcfSSatish Balay 202d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_SMF(Mat mat, Vec v) 203d71ae5a4SJacob Faibussowitsch { 204a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 205a7e14dcfSSatish Balay 206a7e14dcfSSatish Balay PetscFunctionBegin; 2079566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 2089566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(ctx->A, v)); 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 210a7e14dcfSSatish Balay } 211a7e14dcfSSatish Balay 212d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetRowMax_SMF(Mat M, Vec D) 213d71ae5a4SJacob Faibussowitsch { 214a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 215a7e14dcfSSatish Balay 216a7e14dcfSSatish Balay PetscFunctionBegin; 2179566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(M, &ctx)); 2189566063dSJacob Faibussowitsch PetscCall(MatGetRowMax(ctx->A, D, NULL)); 2193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 220a7e14dcfSSatish Balay } 221a7e14dcfSSatish Balay 222d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrices_SMF(Mat A, PetscInt n, IS *irow, IS *icol, MatReuse scall, Mat **B) 223d71ae5a4SJacob Faibussowitsch { 224a7e14dcfSSatish Balay PetscInt i; 225a7e14dcfSSatish Balay 226a7e14dcfSSatish Balay PetscFunctionBegin; 22748a46eb9SPierre Jolivet if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscCalloc1(n + 1, B)); 228a7e14dcfSSatish Balay 22948a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(MatCreateSubMatrix_SMF(A, irow[i], icol[i], scall, &(*B)[i])); 2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 231a7e14dcfSSatish Balay } 232a7e14dcfSSatish Balay 233d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrix_SMF(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat) 234d71ae5a4SJacob Faibussowitsch { 235a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 236a7e14dcfSSatish Balay 237a7e14dcfSSatish Balay PetscFunctionBegin; 2389566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 23948a46eb9SPierre Jolivet if (newmat) PetscCall(MatDestroy(&*newmat)); 2409566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrixFree(ctx->A, isrow, iscol, newmat)); 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 242a7e14dcfSSatish Balay } 243a7e14dcfSSatish Balay 244d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetRow_SMF(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt **cols, const PetscScalar **vals) 245d71ae5a4SJacob Faibussowitsch { 246a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 247a7e14dcfSSatish Balay 248a7e14dcfSSatish Balay PetscFunctionBegin; 2499566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 2509566063dSJacob Faibussowitsch PetscCall(MatGetRow(ctx->A, row, ncols, cols, vals)); 2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 252a7e14dcfSSatish Balay } 253a7e14dcfSSatish Balay 254d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRestoreRow_SMF(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt **cols, const PetscScalar **vals) 255d71ae5a4SJacob Faibussowitsch { 256a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 257a7e14dcfSSatish Balay 258a7e14dcfSSatish Balay PetscFunctionBegin; 2599566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 2609566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(ctx->A, row, ncols, cols, vals)); 2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 262a7e14dcfSSatish Balay } 263a7e14dcfSSatish Balay 264d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnVector_SMF(Mat mat, Vec Y, PetscInt col) 265d71ae5a4SJacob Faibussowitsch { 266a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 267a7e14dcfSSatish Balay 268a7e14dcfSSatish Balay PetscFunctionBegin; 2699566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 2709566063dSJacob Faibussowitsch PetscCall(MatGetColumnVector(ctx->A, Y, col)); 2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 272a7e14dcfSSatish Balay } 273a7e14dcfSSatish Balay 274d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNorm_SMF(Mat mat, NormType type, PetscReal *norm) 275d71ae5a4SJacob Faibussowitsch { 276a7e14dcfSSatish Balay MatSubMatFreeCtx ctx; 277a7e14dcfSSatish Balay 278a7e14dcfSSatish Balay PetscFunctionBegin; 2799566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 280a7e14dcfSSatish Balay if (type == NORM_FROBENIUS) { 281a7e14dcfSSatish Balay *norm = 1.0; 282a7e14dcfSSatish Balay } else if (type == NORM_1 || type == NORM_INFINITY) { 283a7e14dcfSSatish Balay *norm = 1.0; 28487f595a5SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No two norm"); 2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 286a7e14dcfSSatish Balay } 287