xref: /petsc/src/tao/matrix/submatfree.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1ba92ff59SBarry Smith #include <petsctao.h>                     /*I "petsctao.h" I*/
2aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h> /*I "submatfree.h" I*/
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay /*@C
5a7e14dcfSSatish Balay   MatCreateSubMatrixFree - Creates a reduced matrix by masking a
6a7e14dcfSSatish Balay   full matrix.
7a7e14dcfSSatish Balay 
8a7e14dcfSSatish Balay    Collective on matrix
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 
15a7e14dcfSSatish Balay    Output Parameters:
16a7e14dcfSSatish Balay .  J - New matrix
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay    Notes:
190c0fd78eSBarry Smith    The caller is responsible for destroying the input objects after matrix J has been destroyed.
20a7e14dcfSSatish Balay 
21a7e14dcfSSatish Balay    Level: developer
22a7e14dcfSSatish Balay 
23db781477SPatrick Sanan .seealso: `MatCreate()`
24a7e14dcfSSatish Balay @*/
259371c9d4SSatish Balay PetscErrorCode MatCreateSubMatrixFree(Mat mat, IS Rows, IS Cols, Mat *J) {
26367daffbSBarry Smith   MPI_Comm         comm = PetscObjectComm((PetscObject)mat);
27a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
28a7e14dcfSSatish Balay   PetscInt         mloc, nloc, m, n;
29a7e14dcfSSatish Balay 
3087f595a5SBarry Smith   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
32a7e14dcfSSatish Balay   ctx->A = mat;
339566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &m, &n));
349566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mat, &mloc, &nloc));
359566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(mat, NULL, &ctx->VC));
36a7e14dcfSSatish Balay   ctx->VR = ctx->VC;
379566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
38a7e14dcfSSatish Balay 
39a7e14dcfSSatish Balay   ctx->Rows = Rows;
40a7e14dcfSSatish Balay   ctx->Cols = Cols;
419566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Rows));
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Cols));
439566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(comm, mloc, nloc, m, n, ctx, J));
449566063dSJacob Faibussowitsch   PetscCall(MatShellSetManageScalingShifts(*J));
459566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))MatMult_SMF));
469566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))MatDestroy_SMF));
479566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_VIEW, (void (*)(void))MatView_SMF));
489566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_SMF));
499566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DIAGONAL_SET, (void (*)(void))MatDiagonalSet_SMF));
509566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_SHIFT, (void (*)(void))MatShift_SMF));
519566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_EQUAL, (void (*)(void))MatEqual_SMF));
529566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_SCALE, (void (*)(void))MatScale_SMF));
539566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_TRANSPOSE, (void (*)(void))MatTranspose_SMF));
549566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_SMF));
559566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRICES, (void (*)(void))MatCreateSubMatrices_SMF));
569566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_NORM, (void (*)(void))MatNorm_SMF));
579566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_SMF));
589566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRIX, (void (*)(void))MatCreateSubMatrix_SMF));
599566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_GET_ROW_MAX, (void (*)(void))MatDuplicate_SMF));
60a7e14dcfSSatish Balay 
619566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat, (PetscObject)(*J)));
62a7e14dcfSSatish Balay   PetscFunctionReturn(0);
63a7e14dcfSSatish Balay }
64a7e14dcfSSatish Balay 
659371c9d4SSatish Balay PetscErrorCode MatSMFResetRowColumn(Mat mat, IS Rows, IS Cols) {
66a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
67a7e14dcfSSatish Balay 
68a7e14dcfSSatish Balay   PetscFunctionBegin;
699566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
709566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx->Rows));
719566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx->Cols));
729566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Rows));
739566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Cols));
74a7e14dcfSSatish Balay   ctx->Cols = Cols;
75a7e14dcfSSatish Balay   ctx->Rows = Rows;
76a7e14dcfSSatish Balay   PetscFunctionReturn(0);
77a7e14dcfSSatish Balay }
78a7e14dcfSSatish Balay 
799371c9d4SSatish Balay PetscErrorCode MatMult_SMF(Mat mat, Vec a, Vec y) {
80a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
81a7e14dcfSSatish Balay 
82a7e14dcfSSatish Balay   PetscFunctionBegin;
839566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
849566063dSJacob Faibussowitsch   PetscCall(VecCopy(a, ctx->VR));
859566063dSJacob Faibussowitsch   PetscCall(VecISSet(ctx->VR, ctx->Cols, 0.0));
869566063dSJacob Faibussowitsch   PetscCall(MatMult(ctx->A, ctx->VR, y));
879566063dSJacob Faibussowitsch   PetscCall(VecISSet(y, ctx->Rows, 0.0));
88a7e14dcfSSatish Balay   PetscFunctionReturn(0);
89a7e14dcfSSatish Balay }
90a7e14dcfSSatish Balay 
919371c9d4SSatish Balay PetscErrorCode MatMultTranspose_SMF(Mat mat, Vec a, Vec y) {
92a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
969566063dSJacob Faibussowitsch   PetscCall(VecCopy(a, ctx->VC));
979566063dSJacob Faibussowitsch   PetscCall(VecISSet(ctx->VC, ctx->Rows, 0.0));
989566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(ctx->A, ctx->VC, y));
999566063dSJacob Faibussowitsch   PetscCall(VecISSet(y, ctx->Cols, 0.0));
100a7e14dcfSSatish Balay   PetscFunctionReturn(0);
101a7e14dcfSSatish Balay }
102a7e14dcfSSatish Balay 
1039371c9d4SSatish Balay PetscErrorCode MatDiagonalSet_SMF(Mat M, Vec D, InsertMode is) {
104a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
105a7e14dcfSSatish Balay 
106a7e14dcfSSatish Balay   PetscFunctionBegin;
1079566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M, &ctx));
1089566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(ctx->A, D, is));
109a7e14dcfSSatish Balay   PetscFunctionReturn(0);
110a7e14dcfSSatish Balay }
111a7e14dcfSSatish Balay 
1129371c9d4SSatish Balay PetscErrorCode MatDestroy_SMF(Mat mat) {
113a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
114a7e14dcfSSatish Balay 
115a7e14dcfSSatish Balay   PetscFunctionBegin;
1169566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->A));
1189566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx->Rows));
1199566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx->Cols));
1209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->VC));
1219566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
122a7e14dcfSSatish Balay   PetscFunctionReturn(0);
123a7e14dcfSSatish Balay }
124a7e14dcfSSatish Balay 
1259371c9d4SSatish Balay PetscErrorCode MatView_SMF(Mat mat, PetscViewer viewer) {
126a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
127a7e14dcfSSatish Balay 
128a7e14dcfSSatish Balay   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1309566063dSJacob Faibussowitsch   PetscCall(MatView(ctx->A, viewer));
131a7e14dcfSSatish Balay   PetscFunctionReturn(0);
132a7e14dcfSSatish Balay }
133a7e14dcfSSatish Balay 
1349371c9d4SSatish Balay PetscErrorCode MatShift_SMF(Mat Y, PetscReal a) {
135a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
136a7e14dcfSSatish Balay 
137a7e14dcfSSatish Balay   PetscFunctionBegin;
1389566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(Y, &ctx));
1399566063dSJacob Faibussowitsch   PetscCall(MatShift(ctx->A, a));
140a7e14dcfSSatish Balay   PetscFunctionReturn(0);
141a7e14dcfSSatish Balay }
142a7e14dcfSSatish Balay 
1439371c9d4SSatish Balay PetscErrorCode MatDuplicate_SMF(Mat mat, MatDuplicateOption op, Mat *M) {
144a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
145a7e14dcfSSatish Balay 
146a7e14dcfSSatish Balay   PetscFunctionBegin;
1479566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1489566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrixFree(ctx->A, ctx->Rows, ctx->Cols, M));
149a7e14dcfSSatish Balay   PetscFunctionReturn(0);
150a7e14dcfSSatish Balay }
151a7e14dcfSSatish Balay 
1529371c9d4SSatish Balay PetscErrorCode MatEqual_SMF(Mat A, Mat B, PetscBool *flg) {
153a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx1, ctx2;
154a7e14dcfSSatish Balay   PetscBool        flg1, flg2, flg3;
155a7e14dcfSSatish Balay 
156a7e14dcfSSatish Balay   PetscFunctionBegin;
1579566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx1));
1589566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(B, &ctx2));
1599566063dSJacob Faibussowitsch   PetscCall(ISEqual(ctx1->Rows, ctx2->Rows, &flg2));
1609566063dSJacob Faibussowitsch   PetscCall(ISEqual(ctx1->Cols, ctx2->Cols, &flg3));
161a7e14dcfSSatish Balay   if (flg2 == PETSC_FALSE || flg3 == PETSC_FALSE) {
162a7e14dcfSSatish Balay     *flg = PETSC_FALSE;
163a7e14dcfSSatish Balay   } else {
1649566063dSJacob Faibussowitsch     PetscCall(MatEqual(ctx1->A, ctx2->A, &flg1));
1659371c9d4SSatish Balay     if (flg1 == PETSC_FALSE) {
1669371c9d4SSatish Balay       *flg = PETSC_FALSE;
1679371c9d4SSatish Balay     } else {
1689371c9d4SSatish Balay       *flg = PETSC_TRUE;
1699371c9d4SSatish Balay     }
170a7e14dcfSSatish Balay   }
171a7e14dcfSSatish Balay   PetscFunctionReturn(0);
172a7e14dcfSSatish Balay }
173a7e14dcfSSatish Balay 
1749371c9d4SSatish Balay PetscErrorCode MatScale_SMF(Mat mat, PetscReal a) {
175a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
176a7e14dcfSSatish Balay 
177a7e14dcfSSatish Balay   PetscFunctionBegin;
1789566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1799566063dSJacob Faibussowitsch   PetscCall(MatScale(ctx->A, a));
180a7e14dcfSSatish Balay   PetscFunctionReturn(0);
181a7e14dcfSSatish Balay }
182a7e14dcfSSatish Balay 
1839371c9d4SSatish Balay PetscErrorCode MatTranspose_SMF(Mat mat, Mat *B) {
184a7e14dcfSSatish Balay   PetscFunctionBegin;
185a7e14dcfSSatish Balay   PetscFunctionReturn(1);
186a7e14dcfSSatish Balay }
187a7e14dcfSSatish Balay 
1889371c9d4SSatish Balay PetscErrorCode MatGetDiagonal_SMF(Mat mat, Vec v) {
189a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
190a7e14dcfSSatish Balay 
191a7e14dcfSSatish Balay   PetscFunctionBegin;
1929566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1939566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(ctx->A, v));
194a7e14dcfSSatish Balay   PetscFunctionReturn(0);
195a7e14dcfSSatish Balay }
196a7e14dcfSSatish Balay 
1979371c9d4SSatish Balay PetscErrorCode MatGetRowMax_SMF(Mat M, Vec D) {
198a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay   PetscFunctionBegin;
2019566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M, &ctx));
2029566063dSJacob Faibussowitsch   PetscCall(MatGetRowMax(ctx->A, D, NULL));
203a7e14dcfSSatish Balay   PetscFunctionReturn(0);
204a7e14dcfSSatish Balay }
205a7e14dcfSSatish Balay 
2069371c9d4SSatish Balay PetscErrorCode MatCreateSubMatrices_SMF(Mat A, PetscInt n, IS *irow, IS *icol, MatReuse scall, Mat **B) {
207a7e14dcfSSatish Balay   PetscInt i;
208a7e14dcfSSatish Balay 
209a7e14dcfSSatish Balay   PetscFunctionBegin;
210*48a46eb9SPierre Jolivet   if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscCalloc1(n + 1, B));
211a7e14dcfSSatish Balay 
212*48a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(MatCreateSubMatrix_SMF(A, irow[i], icol[i], scall, &(*B)[i]));
213a7e14dcfSSatish Balay   PetscFunctionReturn(0);
214a7e14dcfSSatish Balay }
215a7e14dcfSSatish Balay 
2169371c9d4SSatish Balay PetscErrorCode MatCreateSubMatrix_SMF(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat) {
217a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
218a7e14dcfSSatish Balay 
219a7e14dcfSSatish Balay   PetscFunctionBegin;
2209566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
221*48a46eb9SPierre Jolivet   if (newmat) PetscCall(MatDestroy(&*newmat));
2229566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrixFree(ctx->A, isrow, iscol, newmat));
223a7e14dcfSSatish Balay   PetscFunctionReturn(0);
224a7e14dcfSSatish Balay }
225a7e14dcfSSatish Balay 
2269371c9d4SSatish Balay PetscErrorCode MatGetRow_SMF(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt **cols, const PetscScalar **vals) {
227a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
228a7e14dcfSSatish Balay 
229a7e14dcfSSatish Balay   PetscFunctionBegin;
2309566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2319566063dSJacob Faibussowitsch   PetscCall(MatGetRow(ctx->A, row, ncols, cols, vals));
232a7e14dcfSSatish Balay   PetscFunctionReturn(0);
233a7e14dcfSSatish Balay }
234a7e14dcfSSatish Balay 
2359371c9d4SSatish Balay PetscErrorCode MatRestoreRow_SMF(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt **cols, const PetscScalar **vals) {
236a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay   PetscFunctionBegin;
2399566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2409566063dSJacob Faibussowitsch   PetscCall(MatRestoreRow(ctx->A, row, ncols, cols, vals));
241a7e14dcfSSatish Balay   PetscFunctionReturn(0);
242a7e14dcfSSatish Balay }
243a7e14dcfSSatish Balay 
2449371c9d4SSatish Balay PetscErrorCode MatGetColumnVector_SMF(Mat mat, Vec Y, PetscInt col) {
245a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
246a7e14dcfSSatish Balay 
247a7e14dcfSSatish Balay   PetscFunctionBegin;
2489566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2499566063dSJacob Faibussowitsch   PetscCall(MatGetColumnVector(ctx->A, Y, col));
250a7e14dcfSSatish Balay   PetscFunctionReturn(0);
251a7e14dcfSSatish Balay }
252a7e14dcfSSatish Balay 
2539371c9d4SSatish Balay PetscErrorCode MatNorm_SMF(Mat mat, NormType type, PetscReal *norm) {
254a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
255a7e14dcfSSatish Balay 
256a7e14dcfSSatish Balay   PetscFunctionBegin;
2579566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
258a7e14dcfSSatish Balay   if (type == NORM_FROBENIUS) {
259a7e14dcfSSatish Balay     *norm = 1.0;
260a7e14dcfSSatish Balay   } else if (type == NORM_1 || type == NORM_INFINITY) {
261a7e14dcfSSatish Balay     *norm = 1.0;
26287f595a5SBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No two norm");
263a7e14dcfSSatish Balay   PetscFunctionReturn(0);
264a7e14dcfSSatish Balay }
265