xref: /petsc/src/mat/impls/normal/normm.c (revision 27d8b95c4f79a5cfa41eb48445a6461468fd80f6)
1*27d8b95cSPierre Jolivet #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/
2c8a8475eSBarry Smith 
3c8a8475eSBarry Smith typedef struct {
47cfd0b8cSBarry Smith   Mat A;
5a48507d3SJose E. Roman   Mat D; /* local submatrix for diagonal part */
6*27d8b95cSPierre Jolivet   Vec w;
7c8a8475eSBarry Smith } Mat_Normal;
8c8a8475eSBarry Smith 
966976f2fSJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Normal(Mat A, PetscInt is_max, IS is[], PetscInt ov)
10d71ae5a4SJacob Faibussowitsch {
11*27d8b95cSPierre Jolivet   Mat_Normal *a;
12fa80d070SPierre Jolivet   Mat         pattern;
13fa80d070SPierre Jolivet 
14fa80d070SPierre Jolivet   PetscFunctionBegin;
1508401ef6SPierre Jolivet   PetscCheck(ov >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified");
16*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
179566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(a->A, a->A, NULL, &pattern));
189566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(pattern, MATPRODUCT_AtB));
199566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(pattern));
209566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(pattern));
219566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(pattern, is_max, is, ov));
229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pattern));
233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24fa80d070SPierre Jolivet }
25fa80d070SPierre Jolivet 
2666976f2fSJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Normal(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
27d71ae5a4SJacob Faibussowitsch {
28*27d8b95cSPierre Jolivet   Mat_Normal *a;
29*27d8b95cSPierre Jolivet   Mat         B, *suba;
30fa80d070SPierre Jolivet   IS         *row;
31fa80d070SPierre Jolivet   PetscInt    M;
32fa80d070SPierre Jolivet 
33fa80d070SPierre Jolivet   PetscFunctionBegin;
34*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)mat->data)->zrows && !((Mat_Shell *)mat->data)->zcols, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatZeroRows()/MatZeroRowsColumns() after the SubMatrices creation
35*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)mat->data)->axpy, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatAXPY() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatAXPY() after the SubMatrices creation
36*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)mat->data)->left && !((Mat_Shell *)mat->data)->right, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalScale() after the SubMatrices creation with a SubVector
37*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)mat->data)->dshift, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalSet() after the SubMatrices creation with a SubVector
38*27d8b95cSPierre Jolivet   PetscCheck(irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented");
39*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(mat, &a));
40*27d8b95cSPierre Jolivet   B = a->A;
4148a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
429566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &M, NULL));
439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &row));
449566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0]));
459566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(row[0]));
46fa80d070SPierre Jolivet   for (M = 1; M < n; ++M) row[M] = row[0];
479566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba));
48fa80d070SPierre Jolivet   for (M = 0; M < n; ++M) {
499566063dSJacob Faibussowitsch     PetscCall(MatCreateNormal(suba[M], *submat + M));
50*27d8b95cSPierre Jolivet     ((Mat_Shell *)(*submat)[M]->data)->vscale = ((Mat_Shell *)mat->data)->vscale;
51*27d8b95cSPierre Jolivet     ((Mat_Shell *)(*submat)[M]->data)->vshift = ((Mat_Shell *)mat->data)->vshift;
52fa80d070SPierre Jolivet   }
539566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&row[0]));
549566063dSJacob Faibussowitsch   PetscCall(PetscFree(row));
559566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(n, &suba));
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57fa80d070SPierre Jolivet }
58fa80d070SPierre Jolivet 
5966976f2fSJacob Faibussowitsch static PetscErrorCode MatPermute_Normal(Mat A, IS rowp, IS colp, Mat *B)
60d71ae5a4SJacob Faibussowitsch {
61*27d8b95cSPierre Jolivet   Mat_Normal *a;
62*27d8b95cSPierre Jolivet   Mat         C, Aa;
63fa80d070SPierre Jolivet   IS          row;
64fa80d070SPierre Jolivet 
65fa80d070SPierre Jolivet   PetscFunctionBegin;
66*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)A->data)->zrows && !((Mat_Shell *)A->data)->zcols, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatPermute() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatZeroRows()/MatZeroRowsColumns() after the permutation with a permuted zrows and zcols
67*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)A->data)->axpy, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatPermute() if MatAXPY() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatAXPY() after the permutation with a permuted axpy
68*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)A->data)->left && !((Mat_Shell *)A->data)->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatPermute() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalScale() after the permutation with a permuted left and right
69*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)A->data)->dshift, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatPermute() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalSet() after the permutation with a permuted dshift
7008401ef6SPierre Jolivet   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
71*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
72*27d8b95cSPierre Jolivet   Aa = a->A;
739566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row));
749566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(row));
759566063dSJacob Faibussowitsch   PetscCall(MatPermute(Aa, row, colp, &C));
769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&row));
779566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(C, B));
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
79*27d8b95cSPierre Jolivet   ((Mat_Shell *)(*B)->data)->vscale = ((Mat_Shell *)A->data)->vscale;
80*27d8b95cSPierre Jolivet   ((Mat_Shell *)(*B)->data)->vshift = ((Mat_Shell *)A->data)->vshift;
813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82fa80d070SPierre Jolivet }
83fa80d070SPierre Jolivet 
8466976f2fSJacob Faibussowitsch static PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B)
85d71ae5a4SJacob Faibussowitsch {
86*27d8b95cSPierre Jolivet   Mat_Normal *a;
87fa80d070SPierre Jolivet   Mat         C;
88fa80d070SPierre Jolivet 
89fa80d070SPierre Jolivet   PetscFunctionBegin;
90*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
919566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(a->A, op, &C));
929566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(C, B));
939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
94*27d8b95cSPierre Jolivet   if (op == MAT_COPY_VALUES) PetscCall(MatCopy(A, *B, SAME_NONZERO_PATTERN));
953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
96fa80d070SPierre Jolivet }
97fa80d070SPierre Jolivet 
9866976f2fSJacob Faibussowitsch static PetscErrorCode MatCopy_Normal(Mat A, Mat B, MatStructure str)
99d71ae5a4SJacob Faibussowitsch {
100*27d8b95cSPierre Jolivet   Mat_Normal *a, *b;
101fa80d070SPierre Jolivet 
102fa80d070SPierre Jolivet   PetscFunctionBegin;
103*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
104*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(B, &b));
1059566063dSJacob Faibussowitsch   PetscCall(MatCopy(a->A, b->A, str));
1063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
107fa80d070SPierre Jolivet }
108fa80d070SPierre Jolivet 
10966976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_Normal(Mat N, Vec x, Vec y)
110d71ae5a4SJacob Faibussowitsch {
111*27d8b95cSPierre Jolivet   Mat_Normal *Na;
112c8a8475eSBarry Smith 
113c8a8475eSBarry Smith   PetscFunctionBegin;
114*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
115*27d8b95cSPierre Jolivet   PetscCall(MatMult(Na->A, x, Na->w));
1169566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(Na->A, Na->w, y));
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118db090513SMatthew Knepley }
119db090513SMatthew Knepley 
12066976f2fSJacob Faibussowitsch static PetscErrorCode MatDestroy_Normal(Mat N)
121d71ae5a4SJacob Faibussowitsch {
122*27d8b95cSPierre Jolivet   Mat_Normal *Na;
123c8a8475eSBarry Smith 
124c8a8475eSBarry Smith   PetscFunctionBegin;
125*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
1269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
127a48507d3SJose E. Roman   PetscCall(MatDestroy(&Na->D));
1289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->w));
129*27d8b95cSPierre Jolivet   PetscCall(PetscFree(Na));
1309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMat_C", NULL));
1319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_seqaij_C", NULL));
1329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_mpiaij_C", NULL));
1332f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
1342f0a5c5aSJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_hypre_C", NULL));
1352f0a5c5aSJose E. Roman #endif
1369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_seqdense_C", NULL));
1379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_mpidense_C", NULL));
1389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_dense_C", NULL));
139*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatShellSetContext_C", NULL));
1403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
141c8a8475eSBarry Smith }
142c8a8475eSBarry Smith 
14317a6fd94SBarry Smith /*
14417a6fd94SBarry Smith       Slow, nonscalable version
14517a6fd94SBarry Smith */
14666976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Normal(Mat N, Vec v)
147d71ae5a4SJacob Faibussowitsch {
148*27d8b95cSPierre Jolivet   Mat_Normal        *Na;
149*27d8b95cSPierre Jolivet   Mat                A;
150b24ad042SBarry Smith   PetscInt           i, j, rstart, rend, nnz;
151b24ad042SBarry Smith   const PetscInt    *cols;
15217a6fd94SBarry Smith   PetscScalar       *diag, *work, *values;
153b3cc6726SBarry Smith   const PetscScalar *mvalues;
15417a6fd94SBarry Smith 
15517a6fd94SBarry Smith   PetscFunctionBegin;
156*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
157*27d8b95cSPierre Jolivet   A = Na->A;
1589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work));
1599566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(work, A->cmap->N));
1609566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
16117a6fd94SBarry Smith   for (i = rstart; i < rend; i++) {
1629566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues));
163ad540459SPierre Jolivet     for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * mvalues[j];
1649566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues));
16517a6fd94SBarry Smith   }
1661c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N)));
167d0f46423SBarry Smith   rstart = N->cmap->rstart;
168d0f46423SBarry Smith   rend   = N->cmap->rend;
1699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &values));
1709566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart));
1719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &values));
1729566063dSJacob Faibussowitsch   PetscCall(PetscFree2(diag, work));
1733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17417a6fd94SBarry Smith }
175c8a8475eSBarry Smith 
17666976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_Normal(Mat N, Mat *D)
177d71ae5a4SJacob Faibussowitsch {
178*27d8b95cSPierre Jolivet   Mat_Normal *Na;
179*27d8b95cSPierre Jolivet   Mat         M, A;
180a48507d3SJose E. Roman 
181a48507d3SJose E. Roman   PetscFunctionBegin;
182*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)N->data)->zrows && !((Mat_Shell *)N->data)->zcols, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME
183*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)N->data)->axpy, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatAXPY() has been called on the input Mat");                                            // TODO FIXME
184*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)N->data)->left && !((Mat_Shell *)N->data)->right, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME
185*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)N->data)->dshift, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatDiagonalSet() has been called on the input Mat");                                   // TODO FIXME
186*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
187*27d8b95cSPierre Jolivet   A = Na->A;
188a48507d3SJose E. Roman   PetscCall(MatGetDiagonalBlock(A, &M));
189a48507d3SJose E. Roman   PetscCall(MatCreateNormal(M, &Na->D));
190a48507d3SJose E. Roman   *D = Na->D;
1913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
192a48507d3SJose E. Roman }
193a48507d3SJose E. Roman 
19466976f2fSJacob Faibussowitsch static PetscErrorCode MatNormalGetMat_Normal(Mat A, Mat *M)
195d71ae5a4SJacob Faibussowitsch {
196*27d8b95cSPierre Jolivet   Mat_Normal *Aa;
197fa80d070SPierre Jolivet 
198fa80d070SPierre Jolivet   PetscFunctionBegin;
199*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &Aa));
200fa80d070SPierre Jolivet   *M = Aa->A;
2013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
202fa80d070SPierre Jolivet }
203fa80d070SPierre Jolivet 
204fa80d070SPierre Jolivet /*@
20511a5261eSBarry Smith   MatNormalGetMat - Gets the `Mat` object stored inside a `MATNORMAL`
206fa80d070SPierre Jolivet 
2072ef1f0ffSBarry Smith   Logically Collective
208fa80d070SPierre Jolivet 
209fa80d070SPierre Jolivet   Input Parameter:
21011a5261eSBarry Smith . A - the `MATNORMAL` matrix
211fa80d070SPierre Jolivet 
212fa80d070SPierre Jolivet   Output Parameter:
2132ef1f0ffSBarry Smith . M - the matrix object stored inside `A`
214fa80d070SPierre Jolivet 
215fa80d070SPierre Jolivet   Level: intermediate
216fa80d070SPierre Jolivet 
2171cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatCreateNormal()`
218fa80d070SPierre Jolivet @*/
219d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalGetMat(Mat A, Mat *M)
220d71ae5a4SJacob Faibussowitsch {
221fa80d070SPierre Jolivet   PetscFunctionBegin;
222fa80d070SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
223fa80d070SPierre Jolivet   PetscValidType(A, 1);
2244f572ea9SToby Isaac   PetscAssertPointer(M, 2);
225cac4c232SBarry Smith   PetscUseMethod(A, "MatNormalGetMat_C", (Mat, Mat *), (A, M));
2263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
227fa80d070SPierre Jolivet }
228fa80d070SPierre Jolivet 
22966976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Normal_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
230d71ae5a4SJacob Faibussowitsch {
231*27d8b95cSPierre Jolivet   Mat_Normal *Aa;
232fa80d070SPierre Jolivet   Mat         B;
233fa80d070SPierre Jolivet   PetscInt    m, n, M, N;
234fa80d070SPierre Jolivet 
235fa80d070SPierre Jolivet   PetscFunctionBegin;
236*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &Aa));
2379566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
2389566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
239fa80d070SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
240fa80d070SPierre Jolivet     B = *newmat;
2419566063dSJacob Faibussowitsch     PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B));
242fa80d070SPierre Jolivet   } else {
2439566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B));
2449566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
2459566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(B));
2469566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(B));
2479566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE));
248fa80d070SPierre Jolivet   }
2499566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(B));
250fa80d070SPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
2519566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
252fa80d070SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
2539566063dSJacob Faibussowitsch   PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat));
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
255fa80d070SPierre Jolivet }
256fa80d070SPierre Jolivet 
2572f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
25866976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Normal_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
2592f0a5c5aSJose E. Roman {
2602f0a5c5aSJose E. Roman   PetscFunctionBegin;
2612f0a5c5aSJose E. Roman   if (reuse == MAT_INITIAL_MATRIX) {
2622f0a5c5aSJose E. Roman     PetscCall(MatConvert(A, MATAIJ, reuse, B));
2632f0a5c5aSJose E. Roman     PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B));
2642f0a5c5aSJose E. Roman   } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2662f0a5c5aSJose E. Roman }
2672f0a5c5aSJose E. Roman #endif
2682f0a5c5aSJose E. Roman 
269fa80d070SPierre Jolivet typedef struct {
270fa80d070SPierre Jolivet   Mat work[2];
271fa80d070SPierre Jolivet } Normal_Dense;
272fa80d070SPierre Jolivet 
27366976f2fSJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Normal_Dense(Mat C)
274d71ae5a4SJacob Faibussowitsch {
275fa80d070SPierre Jolivet   Mat           A, B;
276fa80d070SPierre Jolivet   Normal_Dense *contents;
277fa80d070SPierre Jolivet   Mat_Normal   *a;
278fa80d070SPierre Jolivet   PetscScalar  *array;
279fa80d070SPierre Jolivet 
280fa80d070SPierre Jolivet   PetscFunctionBegin;
2818f83a4d9SJacob Faibussowitsch   MatCheckProduct(C, 1);
282fa80d070SPierre Jolivet   A = C->product->A;
283fa80d070SPierre Jolivet   B = C->product->B;
284*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
285fa80d070SPierre Jolivet   contents = (Normal_Dense *)C->product->data;
28628b400f6SJacob Faibussowitsch   PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
287*27d8b95cSPierre Jolivet   if (((Mat_Shell *)A->data)->right) {
2889566063dSJacob Faibussowitsch     PetscCall(MatCopy(B, C, SAME_NONZERO_PATTERN));
289*27d8b95cSPierre Jolivet     PetscCall(MatDiagonalScale(C, ((Mat_Shell *)A->data)->right, NULL));
290fa80d070SPierre Jolivet   }
2919566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(contents->work[0]));
2929566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
2939566063dSJacob Faibussowitsch   PetscCall(MatDensePlaceArray(contents->work[1], array));
2949566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(contents->work[1]));
2959566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
2969566063dSJacob Faibussowitsch   PetscCall(MatDenseResetArray(contents->work[1]));
2979566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
2989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
300*27d8b95cSPierre Jolivet   PetscCall(MatScale(C, ((Mat_Shell *)A->data)->vscale));
3013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
302fa80d070SPierre Jolivet }
303fa80d070SPierre Jolivet 
30466976f2fSJacob Faibussowitsch static PetscErrorCode MatNormal_DenseDestroy(void *ctx)
305d71ae5a4SJacob Faibussowitsch {
306fa80d070SPierre Jolivet   Normal_Dense *contents = (Normal_Dense *)ctx;
307fa80d070SPierre Jolivet 
308fa80d070SPierre Jolivet   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(MatDestroy(contents->work));
3109566063dSJacob Faibussowitsch   PetscCall(MatDestroy(contents->work + 1));
3119566063dSJacob Faibussowitsch   PetscCall(PetscFree(contents));
3123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
313fa80d070SPierre Jolivet }
314fa80d070SPierre Jolivet 
31566976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C)
316d71ae5a4SJacob Faibussowitsch {
317fa80d070SPierre Jolivet   Mat           A, B;
318fa80d070SPierre Jolivet   Normal_Dense *contents = NULL;
319fa80d070SPierre Jolivet   Mat_Normal   *a;
320fa80d070SPierre Jolivet   PetscScalar  *array;
321fa80d070SPierre Jolivet   PetscInt      n, N, m, M;
322fa80d070SPierre Jolivet 
323fa80d070SPierre Jolivet   PetscFunctionBegin;
3248f83a4d9SJacob Faibussowitsch   MatCheckProduct(C, 1);
32528b400f6SJacob Faibussowitsch   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
326fa80d070SPierre Jolivet   A = C->product->A;
327fa80d070SPierre Jolivet   B = C->product->B;
328*27d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
329*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)A->data)->zrows && !((Mat_Shell *)A->data)->zcols, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatProductSymbolic() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME
330*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)A->data)->axpy, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatProductSymbolic() if MatAXPY() has been called on the input Mat");          // TODO FIXME
331*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)A->data)->left, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatProductSymbolic() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME
332*27d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)A->data)->dshift, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatProductSymbolic() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME
3339566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(C, &m, &n));
3349566063dSJacob Faibussowitsch   PetscCall(MatGetSize(C, &M, &N));
335fa80d070SPierre Jolivet   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
3369566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(B, NULL, &n));
3379566063dSJacob Faibussowitsch     PetscCall(MatGetSize(B, NULL, &N));
3389566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &m, NULL));
3399566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &M, NULL));
3409566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, m, n, M, N));
341fa80d070SPierre Jolivet   }
3429566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
3439566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
3449566063dSJacob Faibussowitsch   PetscCall(PetscNew(&contents));
345fa80d070SPierre Jolivet   C->product->data    = contents;
346fa80d070SPierre Jolivet   C->product->destroy = MatNormal_DenseDestroy;
347*27d8b95cSPierre Jolivet   if (((Mat_Shell *)A->data)->right) {
3489566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(a->A, C, NULL, contents->work));
349fa80d070SPierre Jolivet   } else {
3509566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(a->A, B, NULL, contents->work));
351fa80d070SPierre Jolivet   }
3529566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(contents->work[0], MATPRODUCT_AB));
3539566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(contents->work[0]));
3549566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(contents->work[0]));
3559566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(a->A, contents->work[0], NULL, contents->work + 1));
3569566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(contents->work[1], MATPRODUCT_AtB));
3579566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(contents->work[1]));
3589566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(contents->work[1]));
3599566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
3609566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(contents->work[1], array));
3619566063dSJacob Faibussowitsch   PetscCall(MatMPIDenseSetPreallocation(contents->work[1], array));
3629566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
363fa80d070SPierre Jolivet   C->ops->productnumeric = MatProductNumeric_Normal_Dense;
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
365fa80d070SPierre Jolivet }
366fa80d070SPierre Jolivet 
36766976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C)
368d71ae5a4SJacob Faibussowitsch {
369fa80d070SPierre Jolivet   PetscFunctionBegin;
370fa80d070SPierre Jolivet   C->ops->productsymbolic = MatProductSymbolic_Normal_Dense;
3713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
372fa80d070SPierre Jolivet }
373fa80d070SPierre Jolivet 
37466976f2fSJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C)
375d71ae5a4SJacob Faibussowitsch {
376fa80d070SPierre Jolivet   Mat_Product *product = C->product;
377fa80d070SPierre Jolivet 
378fa80d070SPierre Jolivet   PetscFunctionBegin;
37948a46eb9SPierre Jolivet   if (product->type == MATPRODUCT_AB) PetscCall(MatProductSetFromOptions_Normal_Dense_AB(C));
3803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
381fa80d070SPierre Jolivet }
382fa80d070SPierre Jolivet 
3832ef1f0ffSBarry Smith /*MC
3842ef1f0ffSBarry Smith   MATNORMAL - a matrix that behaves like A'*A for `MatMult()` while only containing A
3852ef1f0ffSBarry Smith 
3862ef1f0ffSBarry Smith   Level: intermediate
3872ef1f0ffSBarry Smith 
388*27d8b95cSPierre Jolivet   Developer Notes:
389*27d8b95cSPierre Jolivet   This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
390*27d8b95cSPierre Jolivet 
391*27d8b95cSPierre Jolivet   Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage
392*27d8b95cSPierre Jolivet 
3931cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
3942ef1f0ffSBarry Smith M*/
3952ef1f0ffSBarry Smith 
396c8a8475eSBarry Smith /*@
39711a5261eSBarry Smith   MatCreateNormal - Creates a new `MATNORMAL` matrix object that behaves like A'*A.
398c8a8475eSBarry Smith 
399c3339decSBarry Smith   Collective
400c8a8475eSBarry Smith 
401c8a8475eSBarry Smith   Input Parameter:
402c8a8475eSBarry Smith . A - the (possibly rectangular) matrix
403c8a8475eSBarry Smith 
404c8a8475eSBarry Smith   Output Parameter:
405c8a8475eSBarry Smith . N - the matrix that represents A'*A
406c8a8475eSBarry Smith 
407c30e7cdbSSatish Balay   Level: intermediate
408c30e7cdbSSatish Balay 
40995452b02SPatrick Sanan   Notes:
41095452b02SPatrick Sanan   The product A'*A is NOT actually formed! Rather the new matrix
41111a5261eSBarry Smith   object performs the matrix-vector product, `MatMult()`, by first multiplying by
412c8a8475eSBarry Smith   A and then A'
41311a5261eSBarry Smith 
4141cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
415c8a8475eSBarry Smith @*/
416d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormal(Mat A, Mat *N)
417d71ae5a4SJacob Faibussowitsch {
418c8a8475eSBarry Smith   Mat_Normal *Na;
4199ca0e1b6SStefano Zampini   VecType     vtype;
420c8a8475eSBarry Smith 
421c8a8475eSBarry Smith   PetscFunctionBegin;
4229566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N));
4239566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap));
4249566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap));
425*27d8b95cSPierre Jolivet   PetscCall(MatSetType(*N, MATSHELL));
4264dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&Na));
427*27d8b95cSPierre Jolivet   PetscCall(MatShellSetContext(*N, Na));
4289566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
429c3122656SLisandro Dalcin   Na->A = A;
4309566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, NULL, &Na->w));
4312205254eSKarl Rupp 
432*27d8b95cSPierre Jolivet   PetscCall(MatSetBlockSizes(*N, PetscAbs(A->cmap->bs), PetscAbs(A->rmap->bs)));
433*27d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_DESTROY, (void (*)(void))MatDestroy_Normal));
434*27d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_MULT, (void (*)(void))MatMult_Normal));
435*27d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_Normal));
436*27d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_Normal));
437*27d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Normal));
438*27d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_COPY, (void (*)(void))MatCopy_Normal));
439a48507d3SJose E. Roman   (*N)->ops->getdiagonalblock  = MatGetDiagonalBlock_Normal;
440fa80d070SPierre Jolivet   (*N)->ops->increaseoverlap   = MatIncreaseOverlap_Normal;
441fa80d070SPierre Jolivet   (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal;
442fa80d070SPierre Jolivet   (*N)->ops->permute           = MatPermute_Normal;
4439ca0e1b6SStefano Zampini 
444*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalGetMat_C", MatNormalGetMat_Normal));
445*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_seqaij_C", MatConvert_Normal_AIJ));
446*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_mpiaij_C", MatConvert_Normal_AIJ));
4472f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
448*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_hypre_C", MatConvert_Normal_HYPRE));
4492f0a5c5aSJose E. Roman #endif
450*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_seqdense_C", MatProductSetFromOptions_Normal_Dense));
451*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_mpidense_C", MatProductSetFromOptions_Normal_Dense));
452*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_dense_C", MatProductSetFromOptions_Normal_Dense));
453*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContext_C", MatShellSetContext_Immutable));
454*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
455*27d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
4569566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*N, MAT_SYMMETRIC, PETSC_TRUE));
4579566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A, &vtype));
4589566063dSJacob Faibussowitsch   PetscCall(MatSetVecType(*N, vtype));
4599ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE)
4609566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(*N, A->boundtocpu));
4619ca0e1b6SStefano Zampini #endif
462*27d8b95cSPierre Jolivet   PetscCall(MatSetUp(*N));
463*27d8b95cSPierre Jolivet   PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMAL));
4643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
465c8a8475eSBarry Smith }
466