1 #include <petsctao.h> /*I "petsctao.h" I*/ 2 #include <../src/tao/matrix/submatfree.h> /*I "submatfree.h" I*/ 3 4 /*@C 5 MatCreateSubMatrixFree - Creates a reduced matrix by masking a 6 full matrix. 7 8 Collective on matrix 9 10 Input Parameters: 11 + mat - matrix of arbitrary type 12 . Rows - the rows that will be in the submatrix 13 - Cols - the columns that will be in the submatrix 14 15 Output Parameters: 16 . J - New matrix 17 18 Notes: 19 The caller is responsible for destroying the input objects after matrix J has been destroyed. 20 21 Level: developer 22 23 .seealso: `MatCreate()` 24 @*/ 25 PetscErrorCode MatCreateSubMatrixFree(Mat mat, IS Rows, IS Cols, Mat *J) { 26 MPI_Comm comm = PetscObjectComm((PetscObject)mat); 27 MatSubMatFreeCtx ctx; 28 PetscInt mloc, nloc, m, n; 29 30 PetscFunctionBegin; 31 PetscCall(PetscNew(&ctx)); 32 ctx->A = mat; 33 PetscCall(MatGetSize(mat, &m, &n)); 34 PetscCall(MatGetLocalSize(mat, &mloc, &nloc)); 35 PetscCall(MatCreateVecs(mat, NULL, &ctx->VC)); 36 ctx->VR = ctx->VC; 37 PetscCall(PetscObjectReference((PetscObject)mat)); 38 39 ctx->Rows = Rows; 40 ctx->Cols = Cols; 41 PetscCall(PetscObjectReference((PetscObject)Rows)); 42 PetscCall(PetscObjectReference((PetscObject)Cols)); 43 PetscCall(MatCreateShell(comm, mloc, nloc, m, n, ctx, J)); 44 PetscCall(MatShellSetManageScalingShifts(*J)); 45 PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))MatMult_SMF)); 46 PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))MatDestroy_SMF)); 47 PetscCall(MatShellSetOperation(*J, MATOP_VIEW, (void (*)(void))MatView_SMF)); 48 PetscCall(MatShellSetOperation(*J, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_SMF)); 49 PetscCall(MatShellSetOperation(*J, MATOP_DIAGONAL_SET, (void (*)(void))MatDiagonalSet_SMF)); 50 PetscCall(MatShellSetOperation(*J, MATOP_SHIFT, (void (*)(void))MatShift_SMF)); 51 PetscCall(MatShellSetOperation(*J, MATOP_EQUAL, (void (*)(void))MatEqual_SMF)); 52 PetscCall(MatShellSetOperation(*J, MATOP_SCALE, (void (*)(void))MatScale_SMF)); 53 PetscCall(MatShellSetOperation(*J, MATOP_TRANSPOSE, (void (*)(void))MatTranspose_SMF)); 54 PetscCall(MatShellSetOperation(*J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_SMF)); 55 PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRICES, (void (*)(void))MatCreateSubMatrices_SMF)); 56 PetscCall(MatShellSetOperation(*J, MATOP_NORM, (void (*)(void))MatNorm_SMF)); 57 PetscCall(MatShellSetOperation(*J, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_SMF)); 58 PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRIX, (void (*)(void))MatCreateSubMatrix_SMF)); 59 PetscCall(MatShellSetOperation(*J, MATOP_GET_ROW_MAX, (void (*)(void))MatDuplicate_SMF)); 60 61 PetscFunctionReturn(0); 62 } 63 64 PetscErrorCode MatSMFResetRowColumn(Mat mat, IS Rows, IS Cols) { 65 MatSubMatFreeCtx ctx; 66 67 PetscFunctionBegin; 68 PetscCall(MatShellGetContext(mat, &ctx)); 69 PetscCall(ISDestroy(&ctx->Rows)); 70 PetscCall(ISDestroy(&ctx->Cols)); 71 PetscCall(PetscObjectReference((PetscObject)Rows)); 72 PetscCall(PetscObjectReference((PetscObject)Cols)); 73 ctx->Cols = Cols; 74 ctx->Rows = Rows; 75 PetscFunctionReturn(0); 76 } 77 78 PetscErrorCode MatMult_SMF(Mat mat, Vec a, Vec y) { 79 MatSubMatFreeCtx ctx; 80 81 PetscFunctionBegin; 82 PetscCall(MatShellGetContext(mat, &ctx)); 83 PetscCall(VecCopy(a, ctx->VR)); 84 PetscCall(VecISSet(ctx->VR, ctx->Cols, 0.0)); 85 PetscCall(MatMult(ctx->A, ctx->VR, y)); 86 PetscCall(VecISSet(y, ctx->Rows, 0.0)); 87 PetscFunctionReturn(0); 88 } 89 90 PetscErrorCode MatMultTranspose_SMF(Mat mat, Vec a, Vec y) { 91 MatSubMatFreeCtx ctx; 92 93 PetscFunctionBegin; 94 PetscCall(MatShellGetContext(mat, &ctx)); 95 PetscCall(VecCopy(a, ctx->VC)); 96 PetscCall(VecISSet(ctx->VC, ctx->Rows, 0.0)); 97 PetscCall(MatMultTranspose(ctx->A, ctx->VC, y)); 98 PetscCall(VecISSet(y, ctx->Cols, 0.0)); 99 PetscFunctionReturn(0); 100 } 101 102 PetscErrorCode MatDiagonalSet_SMF(Mat M, Vec D, InsertMode is) { 103 MatSubMatFreeCtx ctx; 104 105 PetscFunctionBegin; 106 PetscCall(MatShellGetContext(M, &ctx)); 107 PetscCall(MatDiagonalSet(ctx->A, D, is)); 108 PetscFunctionReturn(0); 109 } 110 111 PetscErrorCode MatDestroy_SMF(Mat mat) { 112 MatSubMatFreeCtx ctx; 113 114 PetscFunctionBegin; 115 PetscCall(MatShellGetContext(mat, &ctx)); 116 PetscCall(MatDestroy(&ctx->A)); 117 PetscCall(ISDestroy(&ctx->Rows)); 118 PetscCall(ISDestroy(&ctx->Cols)); 119 PetscCall(VecDestroy(&ctx->VC)); 120 PetscCall(PetscFree(ctx)); 121 PetscFunctionReturn(0); 122 } 123 124 PetscErrorCode MatView_SMF(Mat mat, PetscViewer viewer) { 125 MatSubMatFreeCtx ctx; 126 127 PetscFunctionBegin; 128 PetscCall(MatShellGetContext(mat, &ctx)); 129 PetscCall(MatView(ctx->A, viewer)); 130 PetscFunctionReturn(0); 131 } 132 133 PetscErrorCode MatShift_SMF(Mat Y, PetscReal a) { 134 MatSubMatFreeCtx ctx; 135 136 PetscFunctionBegin; 137 PetscCall(MatShellGetContext(Y, &ctx)); 138 PetscCall(MatShift(ctx->A, a)); 139 PetscFunctionReturn(0); 140 } 141 142 PetscErrorCode MatDuplicate_SMF(Mat mat, MatDuplicateOption op, Mat *M) { 143 MatSubMatFreeCtx ctx; 144 145 PetscFunctionBegin; 146 PetscCall(MatShellGetContext(mat, &ctx)); 147 PetscCall(MatCreateSubMatrixFree(ctx->A, ctx->Rows, ctx->Cols, M)); 148 PetscFunctionReturn(0); 149 } 150 151 PetscErrorCode MatEqual_SMF(Mat A, Mat B, PetscBool *flg) { 152 MatSubMatFreeCtx ctx1, ctx2; 153 PetscBool flg1, flg2, flg3; 154 155 PetscFunctionBegin; 156 PetscCall(MatShellGetContext(A, &ctx1)); 157 PetscCall(MatShellGetContext(B, &ctx2)); 158 PetscCall(ISEqual(ctx1->Rows, ctx2->Rows, &flg2)); 159 PetscCall(ISEqual(ctx1->Cols, ctx2->Cols, &flg3)); 160 if (flg2 == PETSC_FALSE || flg3 == PETSC_FALSE) { 161 *flg = PETSC_FALSE; 162 } else { 163 PetscCall(MatEqual(ctx1->A, ctx2->A, &flg1)); 164 if (flg1 == PETSC_FALSE) { 165 *flg = PETSC_FALSE; 166 } else { 167 *flg = PETSC_TRUE; 168 } 169 } 170 PetscFunctionReturn(0); 171 } 172 173 PetscErrorCode MatScale_SMF(Mat mat, PetscReal a) { 174 MatSubMatFreeCtx ctx; 175 176 PetscFunctionBegin; 177 PetscCall(MatShellGetContext(mat, &ctx)); 178 PetscCall(MatScale(ctx->A, a)); 179 PetscFunctionReturn(0); 180 } 181 182 PetscErrorCode MatTranspose_SMF(Mat mat, Mat *B) { 183 PetscFunctionBegin; 184 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "No support for transpose for MatCreateSubMatrixFree() matrix"); 185 PetscFunctionReturn(0); 186 } 187 188 PetscErrorCode MatGetDiagonal_SMF(Mat mat, Vec v) { 189 MatSubMatFreeCtx ctx; 190 191 PetscFunctionBegin; 192 PetscCall(MatShellGetContext(mat, &ctx)); 193 PetscCall(MatGetDiagonal(ctx->A, v)); 194 PetscFunctionReturn(0); 195 } 196 197 PetscErrorCode MatGetRowMax_SMF(Mat M, Vec D) { 198 MatSubMatFreeCtx ctx; 199 200 PetscFunctionBegin; 201 PetscCall(MatShellGetContext(M, &ctx)); 202 PetscCall(MatGetRowMax(ctx->A, D, NULL)); 203 PetscFunctionReturn(0); 204 } 205 206 PetscErrorCode MatCreateSubMatrices_SMF(Mat A, PetscInt n, IS *irow, IS *icol, MatReuse scall, Mat **B) { 207 PetscInt i; 208 209 PetscFunctionBegin; 210 if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscCalloc1(n + 1, B)); 211 212 for (i = 0; i < n; i++) PetscCall(MatCreateSubMatrix_SMF(A, irow[i], icol[i], scall, &(*B)[i])); 213 PetscFunctionReturn(0); 214 } 215 216 PetscErrorCode MatCreateSubMatrix_SMF(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat) { 217 MatSubMatFreeCtx ctx; 218 219 PetscFunctionBegin; 220 PetscCall(MatShellGetContext(mat, &ctx)); 221 if (newmat) PetscCall(MatDestroy(&*newmat)); 222 PetscCall(MatCreateSubMatrixFree(ctx->A, isrow, iscol, newmat)); 223 PetscFunctionReturn(0); 224 } 225 226 PetscErrorCode MatGetRow_SMF(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt **cols, const PetscScalar **vals) { 227 MatSubMatFreeCtx ctx; 228 229 PetscFunctionBegin; 230 PetscCall(MatShellGetContext(mat, &ctx)); 231 PetscCall(MatGetRow(ctx->A, row, ncols, cols, vals)); 232 PetscFunctionReturn(0); 233 } 234 235 PetscErrorCode MatRestoreRow_SMF(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt **cols, const PetscScalar **vals) { 236 MatSubMatFreeCtx ctx; 237 238 PetscFunctionBegin; 239 PetscCall(MatShellGetContext(mat, &ctx)); 240 PetscCall(MatRestoreRow(ctx->A, row, ncols, cols, vals)); 241 PetscFunctionReturn(0); 242 } 243 244 PetscErrorCode MatGetColumnVector_SMF(Mat mat, Vec Y, PetscInt col) { 245 MatSubMatFreeCtx ctx; 246 247 PetscFunctionBegin; 248 PetscCall(MatShellGetContext(mat, &ctx)); 249 PetscCall(MatGetColumnVector(ctx->A, Y, col)); 250 PetscFunctionReturn(0); 251 } 252 253 PetscErrorCode MatNorm_SMF(Mat mat, NormType type, PetscReal *norm) { 254 MatSubMatFreeCtx ctx; 255 256 PetscFunctionBegin; 257 PetscCall(MatShellGetContext(mat, &ctx)); 258 if (type == NORM_FROBENIUS) { 259 *norm = 1.0; 260 } else if (type == NORM_1 || type == NORM_INFINITY) { 261 *norm = 1.0; 262 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No two norm"); 263 PetscFunctionReturn(0); 264 } 265