xref: /petsc/src/tao/matrix/submatfree.c (revision 2fe279fdf3e687a416e4eadb7d3c7a82d60442c6)
1ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/
2954e39ddSJose E. Roman #include <../src/tao/matrix/submatfree.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay /*@C
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 
15*2fe279fdSBarry 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 
23db781477SPatrick Sanan .seealso: `MatCreate()`
24a7e14dcfSSatish Balay @*/
25d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrixFree(Mat mat, IS Rows, IS Cols, Mat *J)
26d71ae5a4SJacob Faibussowitsch {
27367daffbSBarry Smith   MPI_Comm         comm = PetscObjectComm((PetscObject)mat);
28a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
29a7e14dcfSSatish Balay   PetscInt         mloc, nloc, m, n;
30a7e14dcfSSatish Balay 
3187f595a5SBarry Smith   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
33a7e14dcfSSatish Balay   ctx->A = mat;
349566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &m, &n));
359566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mat, &mloc, &nloc));
369566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(mat, NULL, &ctx->VC));
37a7e14dcfSSatish Balay   ctx->VR = ctx->VC;
389566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
39a7e14dcfSSatish Balay 
40a7e14dcfSSatish Balay   ctx->Rows = Rows;
41a7e14dcfSSatish Balay   ctx->Cols = Cols;
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Rows));
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Cols));
449566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(comm, mloc, nloc, m, n, ctx, J));
459566063dSJacob Faibussowitsch   PetscCall(MatShellSetManageScalingShifts(*J));
469566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))MatMult_SMF));
479566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))MatDestroy_SMF));
489566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_VIEW, (void (*)(void))MatView_SMF));
499566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_SMF));
509566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DIAGONAL_SET, (void (*)(void))MatDiagonalSet_SMF));
519566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_SHIFT, (void (*)(void))MatShift_SMF));
529566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_EQUAL, (void (*)(void))MatEqual_SMF));
539566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_SCALE, (void (*)(void))MatScale_SMF));
549566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_TRANSPOSE, (void (*)(void))MatTranspose_SMF));
559566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_SMF));
569566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRICES, (void (*)(void))MatCreateSubMatrices_SMF));
579566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_NORM, (void (*)(void))MatNorm_SMF));
589566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_SMF));
599566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_CREATE_SUBMATRIX, (void (*)(void))MatCreateSubMatrix_SMF));
609566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_GET_ROW_MAX, (void (*)(void))MatDuplicate_SMF));
61a7e14dcfSSatish Balay 
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63a7e14dcfSSatish Balay }
64a7e14dcfSSatish Balay 
65d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSMFResetRowColumn(Mat mat, IS Rows, IS Cols)
66d71ae5a4SJacob Faibussowitsch {
67a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
68a7e14dcfSSatish Balay 
69a7e14dcfSSatish Balay   PetscFunctionBegin;
709566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
719566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx->Rows));
729566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx->Cols));
739566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Rows));
749566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Cols));
75a7e14dcfSSatish Balay   ctx->Cols = Cols;
76a7e14dcfSSatish Balay   ctx->Rows = Rows;
773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78a7e14dcfSSatish Balay }
79a7e14dcfSSatish Balay 
80d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_SMF(Mat mat, Vec a, Vec y)
81d71ae5a4SJacob Faibussowitsch {
82a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
83a7e14dcfSSatish Balay 
84a7e14dcfSSatish Balay   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
869566063dSJacob Faibussowitsch   PetscCall(VecCopy(a, ctx->VR));
879566063dSJacob Faibussowitsch   PetscCall(VecISSet(ctx->VR, ctx->Cols, 0.0));
889566063dSJacob Faibussowitsch   PetscCall(MatMult(ctx->A, ctx->VR, y));
899566063dSJacob Faibussowitsch   PetscCall(VecISSet(y, ctx->Rows, 0.0));
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91a7e14dcfSSatish Balay }
92a7e14dcfSSatish Balay 
93d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTranspose_SMF(Mat mat, Vec a, Vec y)
94d71ae5a4SJacob Faibussowitsch {
95a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
999566063dSJacob Faibussowitsch   PetscCall(VecCopy(a, ctx->VC));
1009566063dSJacob Faibussowitsch   PetscCall(VecISSet(ctx->VC, ctx->Rows, 0.0));
1019566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(ctx->A, ctx->VC, y));
1029566063dSJacob Faibussowitsch   PetscCall(VecISSet(y, ctx->Cols, 0.0));
1033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104a7e14dcfSSatish Balay }
105a7e14dcfSSatish Balay 
106d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet_SMF(Mat M, Vec D, InsertMode is)
107d71ae5a4SJacob Faibussowitsch {
108a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
109a7e14dcfSSatish Balay 
110a7e14dcfSSatish Balay   PetscFunctionBegin;
1119566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M, &ctx));
1129566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(ctx->A, D, is));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114a7e14dcfSSatish Balay }
115a7e14dcfSSatish Balay 
116d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_SMF(Mat mat)
117d71ae5a4SJacob Faibussowitsch {
118a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
119a7e14dcfSSatish Balay 
120a7e14dcfSSatish Balay   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->A));
1239566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx->Rows));
1249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx->Cols));
1259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->VC));
1269566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
1273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
128a7e14dcfSSatish Balay }
129a7e14dcfSSatish Balay 
130d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SMF(Mat mat, PetscViewer viewer)
131d71ae5a4SJacob Faibussowitsch {
132a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
133a7e14dcfSSatish Balay 
134a7e14dcfSSatish Balay   PetscFunctionBegin;
1359566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1369566063dSJacob Faibussowitsch   PetscCall(MatView(ctx->A, viewer));
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138a7e14dcfSSatish Balay }
139a7e14dcfSSatish Balay 
140d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift_SMF(Mat Y, PetscReal a)
141d71ae5a4SJacob Faibussowitsch {
142a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
143a7e14dcfSSatish Balay 
144a7e14dcfSSatish Balay   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(Y, &ctx));
1469566063dSJacob Faibussowitsch   PetscCall(MatShift(ctx->A, a));
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
148a7e14dcfSSatish Balay }
149a7e14dcfSSatish Balay 
150d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_SMF(Mat mat, MatDuplicateOption op, Mat *M)
151d71ae5a4SJacob Faibussowitsch {
152a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
153a7e14dcfSSatish Balay 
154a7e14dcfSSatish Balay   PetscFunctionBegin;
1559566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1569566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrixFree(ctx->A, ctx->Rows, ctx->Cols, M));
1573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
158a7e14dcfSSatish Balay }
159a7e14dcfSSatish Balay 
160d71ae5a4SJacob Faibussowitsch PetscErrorCode MatEqual_SMF(Mat A, Mat B, PetscBool *flg)
161d71ae5a4SJacob Faibussowitsch {
162a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx1, ctx2;
163a7e14dcfSSatish Balay   PetscBool        flg1, flg2, flg3;
164a7e14dcfSSatish Balay 
165a7e14dcfSSatish Balay   PetscFunctionBegin;
1669566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx1));
1679566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(B, &ctx2));
1689566063dSJacob Faibussowitsch   PetscCall(ISEqual(ctx1->Rows, ctx2->Rows, &flg2));
1699566063dSJacob Faibussowitsch   PetscCall(ISEqual(ctx1->Cols, ctx2->Cols, &flg3));
170a7e14dcfSSatish Balay   if (flg2 == PETSC_FALSE || flg3 == PETSC_FALSE) {
171a7e14dcfSSatish Balay     *flg = PETSC_FALSE;
172a7e14dcfSSatish Balay   } else {
1739566063dSJacob Faibussowitsch     PetscCall(MatEqual(ctx1->A, ctx2->A, &flg1));
1749371c9d4SSatish Balay     if (flg1 == PETSC_FALSE) {
1759371c9d4SSatish Balay       *flg = PETSC_FALSE;
1769371c9d4SSatish Balay     } else {
1779371c9d4SSatish Balay       *flg = PETSC_TRUE;
1789371c9d4SSatish Balay     }
179a7e14dcfSSatish Balay   }
1803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
181a7e14dcfSSatish Balay }
182a7e14dcfSSatish Balay 
183d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_SMF(Mat mat, PetscReal a)
184d71ae5a4SJacob Faibussowitsch {
185a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
186a7e14dcfSSatish Balay 
187a7e14dcfSSatish Balay   PetscFunctionBegin;
1889566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
1899566063dSJacob Faibussowitsch   PetscCall(MatScale(ctx->A, a));
1903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
191a7e14dcfSSatish Balay }
192a7e14dcfSSatish Balay 
193d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTranspose_SMF(Mat mat, Mat *B)
194d71ae5a4SJacob Faibussowitsch {
195a7e14dcfSSatish Balay   PetscFunctionBegin;
19611cc89d2SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "No support for transpose for MatCreateSubMatrixFree() matrix");
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
198a7e14dcfSSatish Balay }
199a7e14dcfSSatish Balay 
200d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_SMF(Mat mat, Vec v)
201d71ae5a4SJacob Faibussowitsch {
202a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
203a7e14dcfSSatish Balay 
204a7e14dcfSSatish Balay   PetscFunctionBegin;
2059566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2069566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(ctx->A, v));
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
208a7e14dcfSSatish Balay }
209a7e14dcfSSatish Balay 
210d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetRowMax_SMF(Mat M, Vec D)
211d71ae5a4SJacob Faibussowitsch {
212a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
213a7e14dcfSSatish Balay 
214a7e14dcfSSatish Balay   PetscFunctionBegin;
2159566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(M, &ctx));
2169566063dSJacob Faibussowitsch   PetscCall(MatGetRowMax(ctx->A, D, NULL));
2173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
218a7e14dcfSSatish Balay }
219a7e14dcfSSatish Balay 
220d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrices_SMF(Mat A, PetscInt n, IS *irow, IS *icol, MatReuse scall, Mat **B)
221d71ae5a4SJacob Faibussowitsch {
222a7e14dcfSSatish Balay   PetscInt i;
223a7e14dcfSSatish Balay 
224a7e14dcfSSatish Balay   PetscFunctionBegin;
22548a46eb9SPierre Jolivet   if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscCalloc1(n + 1, B));
226a7e14dcfSSatish Balay 
22748a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(MatCreateSubMatrix_SMF(A, irow[i], icol[i], scall, &(*B)[i]));
2283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
229a7e14dcfSSatish Balay }
230a7e14dcfSSatish Balay 
231d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrix_SMF(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat)
232d71ae5a4SJacob Faibussowitsch {
233a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
234a7e14dcfSSatish Balay 
235a7e14dcfSSatish Balay   PetscFunctionBegin;
2369566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
23748a46eb9SPierre Jolivet   if (newmat) PetscCall(MatDestroy(&*newmat));
2389566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrixFree(ctx->A, isrow, iscol, newmat));
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240a7e14dcfSSatish Balay }
241a7e14dcfSSatish Balay 
242d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetRow_SMF(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt **cols, const PetscScalar **vals)
243d71ae5a4SJacob Faibussowitsch {
244a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
245a7e14dcfSSatish Balay 
246a7e14dcfSSatish Balay   PetscFunctionBegin;
2479566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2489566063dSJacob Faibussowitsch   PetscCall(MatGetRow(ctx->A, row, ncols, cols, vals));
2493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
250a7e14dcfSSatish Balay }
251a7e14dcfSSatish Balay 
252d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRestoreRow_SMF(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt **cols, const PetscScalar **vals)
253d71ae5a4SJacob Faibussowitsch {
254a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
255a7e14dcfSSatish Balay 
256a7e14dcfSSatish Balay   PetscFunctionBegin;
2579566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2589566063dSJacob Faibussowitsch   PetscCall(MatRestoreRow(ctx->A, row, ncols, cols, vals));
2593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
260a7e14dcfSSatish Balay }
261a7e14dcfSSatish Balay 
262d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnVector_SMF(Mat mat, Vec Y, PetscInt col)
263d71ae5a4SJacob Faibussowitsch {
264a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
265a7e14dcfSSatish Balay 
266a7e14dcfSSatish Balay   PetscFunctionBegin;
2679566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
2689566063dSJacob Faibussowitsch   PetscCall(MatGetColumnVector(ctx->A, Y, col));
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
270a7e14dcfSSatish Balay }
271a7e14dcfSSatish Balay 
272d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNorm_SMF(Mat mat, NormType type, PetscReal *norm)
273d71ae5a4SJacob Faibussowitsch {
274a7e14dcfSSatish Balay   MatSubMatFreeCtx ctx;
275a7e14dcfSSatish Balay 
276a7e14dcfSSatish Balay   PetscFunctionBegin;
2779566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
278a7e14dcfSSatish Balay   if (type == NORM_FROBENIUS) {
279a7e14dcfSSatish Balay     *norm = 1.0;
280a7e14dcfSSatish Balay   } else if (type == NORM_1 || type == NORM_INFINITY) {
281a7e14dcfSSatish Balay     *norm = 1.0;
28287f595a5SBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No two norm");
2833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
284a7e14dcfSSatish Balay }
285