xref: /petsc/src/mat/impls/normal/normm.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
1c8a8475eSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
3c8a8475eSBarry Smith 
4c8a8475eSBarry Smith typedef struct {
57cfd0b8cSBarry Smith   Mat         A;
6a48507d3SJose E. Roman   Mat         D; /* local submatrix for diagonal part */
77cfd0b8cSBarry Smith   Vec         w, left, right, leftwork, rightwork;
8b20f02adSBarry Smith   PetscScalar scale;
9c8a8475eSBarry Smith } Mat_Normal;
10c8a8475eSBarry Smith 
11d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_Normal(Mat inA, PetscScalar scale)
12d71ae5a4SJacob Faibussowitsch {
1369ef1043SBarry Smith   Mat_Normal *a = (Mat_Normal *)inA->data;
1469ef1043SBarry Smith 
1569ef1043SBarry Smith   PetscFunctionBegin;
1669ef1043SBarry Smith   a->scale *= scale;
173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1869ef1043SBarry Smith }
1969ef1043SBarry Smith 
20d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalScale_Normal(Mat inA, Vec left, Vec right)
21d71ae5a4SJacob Faibussowitsch {
227cfd0b8cSBarry Smith   Mat_Normal *a = (Mat_Normal *)inA->data;
237cfd0b8cSBarry Smith 
247cfd0b8cSBarry Smith   PetscFunctionBegin;
257cfd0b8cSBarry Smith   if (left) {
267cfd0b8cSBarry Smith     if (!a->left) {
279566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(left, &a->left));
289566063dSJacob Faibussowitsch       PetscCall(VecCopy(left, a->left));
297cfd0b8cSBarry Smith     } else {
309566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->left, left, a->left));
317cfd0b8cSBarry Smith     }
327cfd0b8cSBarry Smith   }
337cfd0b8cSBarry Smith   if (right) {
347cfd0b8cSBarry Smith     if (!a->right) {
359566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(right, &a->right));
369566063dSJacob Faibussowitsch       PetscCall(VecCopy(right, a->right));
377cfd0b8cSBarry Smith     } else {
389566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->right, right, a->right));
397cfd0b8cSBarry Smith     }
407cfd0b8cSBarry Smith   }
413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
427cfd0b8cSBarry Smith }
437cfd0b8cSBarry Smith 
44d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIncreaseOverlap_Normal(Mat A, PetscInt is_max, IS is[], PetscInt ov)
45d71ae5a4SJacob Faibussowitsch {
46fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
47fa80d070SPierre Jolivet   Mat         pattern;
48fa80d070SPierre Jolivet 
49fa80d070SPierre Jolivet   PetscFunctionBegin;
5008401ef6SPierre Jolivet   PetscCheck(ov >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified");
519566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(a->A, a->A, NULL, &pattern));
529566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(pattern, MATPRODUCT_AtB));
539566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(pattern));
549566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(pattern));
559566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(pattern, is_max, is, ov));
569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pattern));
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58fa80d070SPierre Jolivet }
59fa80d070SPierre Jolivet 
60d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrices_Normal(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
61d71ae5a4SJacob Faibussowitsch {
62fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)mat->data;
63fa80d070SPierre Jolivet   Mat         B = a->A, *suba;
64fa80d070SPierre Jolivet   IS         *row;
65fa80d070SPierre Jolivet   PetscInt    M;
66fa80d070SPierre Jolivet 
67fa80d070SPierre Jolivet   PetscFunctionBegin;
68aed4548fSBarry Smith   PetscCheck(!a->left && !a->right && irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented");
6948a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
709566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &M, NULL));
719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &row));
729566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0]));
739566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(row[0]));
74fa80d070SPierre Jolivet   for (M = 1; M < n; ++M) row[M] = row[0];
759566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba));
76fa80d070SPierre Jolivet   for (M = 0; M < n; ++M) {
779566063dSJacob Faibussowitsch     PetscCall(MatCreateNormal(suba[M], *submat + M));
78fa80d070SPierre Jolivet     ((Mat_Normal *)(*submat)[M]->data)->scale = a->scale;
79fa80d070SPierre Jolivet   }
809566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&row[0]));
819566063dSJacob Faibussowitsch   PetscCall(PetscFree(row));
829566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(n, &suba));
833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84fa80d070SPierre Jolivet }
85fa80d070SPierre Jolivet 
86d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPermute_Normal(Mat A, IS rowp, IS colp, Mat *B)
87d71ae5a4SJacob Faibussowitsch {
88fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
89fa80d070SPierre Jolivet   Mat         C, Aa = a->A;
90fa80d070SPierre Jolivet   IS          row;
91fa80d070SPierre Jolivet 
92fa80d070SPierre Jolivet   PetscFunctionBegin;
9308401ef6SPierre Jolivet   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
949566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row));
959566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(row));
969566063dSJacob Faibussowitsch   PetscCall(MatPermute(Aa, row, colp, &C));
979566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&row));
989566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(C, B));
999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
101fa80d070SPierre Jolivet }
102fa80d070SPierre Jolivet 
103d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B)
104d71ae5a4SJacob Faibussowitsch {
105fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
106fa80d070SPierre Jolivet   Mat         C;
107fa80d070SPierre Jolivet 
108fa80d070SPierre Jolivet   PetscFunctionBegin;
10908401ef6SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
1109566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(a->A, op, &C));
1119566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(C, B));
1129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
113fa80d070SPierre Jolivet   if (op == MAT_COPY_VALUES) ((Mat_Normal *)(*B)->data)->scale = a->scale;
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115fa80d070SPierre Jolivet }
116fa80d070SPierre Jolivet 
117d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCopy_Normal(Mat A, Mat B, MatStructure str)
118d71ae5a4SJacob Faibussowitsch {
119fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data, *b = (Mat_Normal *)B->data;
120fa80d070SPierre Jolivet 
121fa80d070SPierre Jolivet   PetscFunctionBegin;
12208401ef6SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
1239566063dSJacob Faibussowitsch   PetscCall(MatCopy(a->A, b->A, str));
124fa80d070SPierre Jolivet   b->scale = a->scale;
1259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->left));
1269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->right));
1279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->leftwork));
1289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->rightwork));
1293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
130fa80d070SPierre Jolivet }
131fa80d070SPierre Jolivet 
132d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Normal(Mat N, Vec x, Vec y)
133d71ae5a4SJacob Faibussowitsch {
134c8a8475eSBarry Smith   Mat_Normal *Na = (Mat_Normal *)N->data;
1357cfd0b8cSBarry Smith   Vec         in;
136c8a8475eSBarry Smith 
137c8a8475eSBarry Smith   PetscFunctionBegin;
1387cfd0b8cSBarry Smith   in = x;
1397cfd0b8cSBarry Smith   if (Na->right) {
14048a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1419566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
1427cfd0b8cSBarry Smith     in = Na->rightwork;
1437cfd0b8cSBarry Smith   }
1449566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1459566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(Na->A, Na->w, y));
1461baa6e33SBarry Smith   if (Na->left) PetscCall(VecPointwiseMult(y, Na->left, y));
1479566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
149c8a8475eSBarry Smith }
150c8a8475eSBarry Smith 
151d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3)
152d71ae5a4SJacob Faibussowitsch {
153db090513SMatthew Knepley   Mat_Normal *Na = (Mat_Normal *)N->data;
1547cfd0b8cSBarry Smith   Vec         in;
155db090513SMatthew Knepley 
156db090513SMatthew Knepley   PetscFunctionBegin;
1577cfd0b8cSBarry Smith   in = v1;
1587cfd0b8cSBarry Smith   if (Na->right) {
15948a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1609566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
1617cfd0b8cSBarry Smith     in = Na->rightwork;
1627cfd0b8cSBarry Smith   }
1639566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1649566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
1657cfd0b8cSBarry Smith   if (Na->left) {
1669566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(Na->A, Na->w, v3));
1679566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->left, v3));
1689566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
1697cfd0b8cSBarry Smith   } else {
1709566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(Na->A, Na->w, v2, v3));
1717cfd0b8cSBarry Smith   }
1723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173b20f02adSBarry Smith }
174b20f02adSBarry Smith 
175d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTranspose_Normal(Mat N, Vec x, Vec y)
176d71ae5a4SJacob Faibussowitsch {
177b20f02adSBarry Smith   Mat_Normal *Na = (Mat_Normal *)N->data;
1787cfd0b8cSBarry Smith   Vec         in;
179b20f02adSBarry Smith 
180b20f02adSBarry Smith   PetscFunctionBegin;
1817cfd0b8cSBarry Smith   in = x;
1827cfd0b8cSBarry Smith   if (Na->left) {
18348a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
1849566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
1857cfd0b8cSBarry Smith     in = Na->leftwork;
1867cfd0b8cSBarry Smith   }
1879566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1889566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(Na->A, Na->w, y));
1891baa6e33SBarry Smith   if (Na->right) PetscCall(VecPointwiseMult(y, Na->right, y));
1909566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
1913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
192b20f02adSBarry Smith }
193b20f02adSBarry Smith 
194d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3)
195d71ae5a4SJacob Faibussowitsch {
196b20f02adSBarry Smith   Mat_Normal *Na = (Mat_Normal *)N->data;
1977cfd0b8cSBarry Smith   Vec         in;
198b20f02adSBarry Smith 
199b20f02adSBarry Smith   PetscFunctionBegin;
2007cfd0b8cSBarry Smith   in = v1;
2017cfd0b8cSBarry Smith   if (Na->left) {
20248a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
2039566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
2047cfd0b8cSBarry Smith     in = Na->leftwork;
2057cfd0b8cSBarry Smith   }
2069566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
2079566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
2087cfd0b8cSBarry Smith   if (Na->right) {
2099566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(Na->A, Na->w, v3));
2109566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->right, v3));
2119566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
2127cfd0b8cSBarry Smith   } else {
2139566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(Na->A, Na->w, v2, v3));
2147cfd0b8cSBarry Smith   }
2153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
216db090513SMatthew Knepley }
217db090513SMatthew Knepley 
218d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_Normal(Mat N)
219d71ae5a4SJacob Faibussowitsch {
220c8a8475eSBarry Smith   Mat_Normal *Na = (Mat_Normal *)N->data;
221c8a8475eSBarry Smith 
222c8a8475eSBarry Smith   PetscFunctionBegin;
2239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
224a48507d3SJose E. Roman   PetscCall(MatDestroy(&Na->D));
2259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->w));
2269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->left));
2279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->right));
2289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->leftwork));
2299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->rightwork));
2309566063dSJacob Faibussowitsch   PetscCall(PetscFree(N->data));
2319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMat_C", NULL));
2329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_seqaij_C", NULL));
2339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_mpiaij_C", NULL));
2342f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
2352f0a5c5aSJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_hypre_C", NULL));
2362f0a5c5aSJose E. Roman #endif
2379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_seqdense_C", NULL));
2389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_mpidense_C", NULL));
2399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_dense_C", NULL));
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
241c8a8475eSBarry Smith }
242c8a8475eSBarry Smith 
24317a6fd94SBarry Smith /*
24417a6fd94SBarry Smith       Slow, nonscalable version
24517a6fd94SBarry Smith */
246d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_Normal(Mat N, Vec v)
247d71ae5a4SJacob Faibussowitsch {
24817a6fd94SBarry Smith   Mat_Normal        *Na = (Mat_Normal *)N->data;
24917a6fd94SBarry Smith   Mat                A  = Na->A;
250b24ad042SBarry Smith   PetscInt           i, j, rstart, rend, nnz;
251b24ad042SBarry Smith   const PetscInt    *cols;
25217a6fd94SBarry Smith   PetscScalar       *diag, *work, *values;
253b3cc6726SBarry Smith   const PetscScalar *mvalues;
25417a6fd94SBarry Smith 
25517a6fd94SBarry Smith   PetscFunctionBegin;
2569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work));
2579566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(work, A->cmap->N));
2589566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
25917a6fd94SBarry Smith   for (i = rstart; i < rend; i++) {
2609566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues));
261ad540459SPierre Jolivet     for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * mvalues[j];
2629566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues));
26317a6fd94SBarry Smith   }
2641c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N)));
265d0f46423SBarry Smith   rstart = N->cmap->rstart;
266d0f46423SBarry Smith   rend   = N->cmap->rend;
2679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &values));
2689566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart));
2699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &values));
2709566063dSJacob Faibussowitsch   PetscCall(PetscFree2(diag, work));
2719566063dSJacob Faibussowitsch   PetscCall(VecScale(v, Na->scale));
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27317a6fd94SBarry Smith }
274c8a8475eSBarry Smith 
275d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonalBlock_Normal(Mat N, Mat *D)
276d71ae5a4SJacob Faibussowitsch {
277a48507d3SJose E. Roman   Mat_Normal *Na = (Mat_Normal *)N->data;
278a48507d3SJose E. Roman   Mat         M, A = Na->A;
279a48507d3SJose E. Roman 
280a48507d3SJose E. Roman   PetscFunctionBegin;
281a48507d3SJose E. Roman   PetscCall(MatGetDiagonalBlock(A, &M));
282a48507d3SJose E. Roman   PetscCall(MatCreateNormal(M, &Na->D));
283a48507d3SJose E. Roman   *D = Na->D;
2843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
285a48507d3SJose E. Roman }
286a48507d3SJose E. Roman 
287d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalGetMat_Normal(Mat A, Mat *M)
288d71ae5a4SJacob Faibussowitsch {
289fa80d070SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal *)A->data;
290fa80d070SPierre Jolivet 
291fa80d070SPierre Jolivet   PetscFunctionBegin;
292fa80d070SPierre Jolivet   *M = Aa->A;
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
294fa80d070SPierre Jolivet }
295fa80d070SPierre Jolivet 
296fa80d070SPierre Jolivet /*@
29711a5261eSBarry Smith       MatNormalGetMat - Gets the `Mat` object stored inside a `MATNORMAL`
298fa80d070SPierre Jolivet 
2992ef1f0ffSBarry Smith    Logically Collective
300fa80d070SPierre Jolivet 
301fa80d070SPierre Jolivet    Input Parameter:
30211a5261eSBarry Smith .   A  - the `MATNORMAL` matrix
303fa80d070SPierre Jolivet 
304fa80d070SPierre Jolivet    Output Parameter:
3052ef1f0ffSBarry Smith .   M - the matrix object stored inside `A`
306fa80d070SPierre Jolivet 
307fa80d070SPierre Jolivet    Level: intermediate
308fa80d070SPierre Jolivet 
309*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatCreateNormal()`
310fa80d070SPierre Jolivet @*/
311d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalGetMat(Mat A, Mat *M)
312d71ae5a4SJacob Faibussowitsch {
313fa80d070SPierre Jolivet   PetscFunctionBegin;
314fa80d070SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
315fa80d070SPierre Jolivet   PetscValidType(A, 1);
316fa80d070SPierre Jolivet   PetscValidPointer(M, 2);
317cac4c232SBarry Smith   PetscUseMethod(A, "MatNormalGetMat_C", (Mat, Mat *), (A, M));
3183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
319fa80d070SPierre Jolivet }
320fa80d070SPierre Jolivet 
321d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_Normal_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
322d71ae5a4SJacob Faibussowitsch {
323fa80d070SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal *)A->data;
324fa80d070SPierre Jolivet   Mat         B;
325fa80d070SPierre Jolivet   PetscInt    m, n, M, N;
326fa80d070SPierre Jolivet 
327fa80d070SPierre Jolivet   PetscFunctionBegin;
3289566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
3299566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
330fa80d070SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
331fa80d070SPierre Jolivet     B = *newmat;
3329566063dSJacob Faibussowitsch     PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B));
333fa80d070SPierre Jolivet   } else {
3349566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B));
3359566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
3369566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(B));
3379566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(B));
3389566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE));
339fa80d070SPierre Jolivet   }
3409566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(B));
341fa80d070SPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
3429566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
343fa80d070SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
3449566063dSJacob Faibussowitsch   PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat));
3453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
346fa80d070SPierre Jolivet }
347fa80d070SPierre Jolivet 
3482f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
3492f0a5c5aSJose E. Roman PetscErrorCode MatConvert_Normal_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
3502f0a5c5aSJose E. Roman {
3512f0a5c5aSJose E. Roman   PetscFunctionBegin;
3522f0a5c5aSJose E. Roman   if (reuse == MAT_INITIAL_MATRIX) {
3532f0a5c5aSJose E. Roman     PetscCall(MatConvert(A, MATAIJ, reuse, B));
3542f0a5c5aSJose E. Roman     PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B));
3552f0a5c5aSJose E. Roman   } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */
3563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3572f0a5c5aSJose E. Roman }
3582f0a5c5aSJose E. Roman #endif
3592f0a5c5aSJose E. Roman 
360fa80d070SPierre Jolivet typedef struct {
361fa80d070SPierre Jolivet   Mat work[2];
362fa80d070SPierre Jolivet } Normal_Dense;
363fa80d070SPierre Jolivet 
364d71ae5a4SJacob Faibussowitsch PetscErrorCode MatProductNumeric_Normal_Dense(Mat C)
365d71ae5a4SJacob Faibussowitsch {
366fa80d070SPierre Jolivet   Mat           A, B;
367fa80d070SPierre Jolivet   Normal_Dense *contents;
368fa80d070SPierre Jolivet   Mat_Normal   *a;
369fa80d070SPierre Jolivet   PetscScalar  *array;
370fa80d070SPierre Jolivet 
371fa80d070SPierre Jolivet   PetscFunctionBegin;
3728f83a4d9SJacob Faibussowitsch   MatCheckProduct(C, 1);
373fa80d070SPierre Jolivet   A        = C->product->A;
374fa80d070SPierre Jolivet   a        = (Mat_Normal *)A->data;
375fa80d070SPierre Jolivet   B        = C->product->B;
376fa80d070SPierre Jolivet   contents = (Normal_Dense *)C->product->data;
37728b400f6SJacob Faibussowitsch   PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
378fa80d070SPierre Jolivet   if (a->right) {
3799566063dSJacob Faibussowitsch     PetscCall(MatCopy(B, C, SAME_NONZERO_PATTERN));
3809566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(C, a->right, NULL));
381fa80d070SPierre Jolivet   }
3829566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(contents->work[0]));
3839566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
3849566063dSJacob Faibussowitsch   PetscCall(MatDensePlaceArray(contents->work[1], array));
3859566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(contents->work[1]));
3869566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
3879566063dSJacob Faibussowitsch   PetscCall(MatDenseResetArray(contents->work[1]));
3889566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
3899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
3909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
3919566063dSJacob Faibussowitsch   PetscCall(MatScale(C, a->scale));
3923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
393fa80d070SPierre Jolivet }
394fa80d070SPierre Jolivet 
395d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormal_DenseDestroy(void *ctx)
396d71ae5a4SJacob Faibussowitsch {
397fa80d070SPierre Jolivet   Normal_Dense *contents = (Normal_Dense *)ctx;
398fa80d070SPierre Jolivet 
399fa80d070SPierre Jolivet   PetscFunctionBegin;
4009566063dSJacob Faibussowitsch   PetscCall(MatDestroy(contents->work));
4019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(contents->work + 1));
4029566063dSJacob Faibussowitsch   PetscCall(PetscFree(contents));
4033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
404fa80d070SPierre Jolivet }
405fa80d070SPierre Jolivet 
406d71ae5a4SJacob Faibussowitsch PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C)
407d71ae5a4SJacob Faibussowitsch {
408fa80d070SPierre Jolivet   Mat           A, B;
409fa80d070SPierre Jolivet   Normal_Dense *contents = NULL;
410fa80d070SPierre Jolivet   Mat_Normal   *a;
411fa80d070SPierre Jolivet   PetscScalar  *array;
412fa80d070SPierre Jolivet   PetscInt      n, N, m, M;
413fa80d070SPierre Jolivet 
414fa80d070SPierre Jolivet   PetscFunctionBegin;
4158f83a4d9SJacob Faibussowitsch   MatCheckProduct(C, 1);
41628b400f6SJacob Faibussowitsch   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
417fa80d070SPierre Jolivet   A = C->product->A;
418fa80d070SPierre Jolivet   a = (Mat_Normal *)A->data;
41928b400f6SJacob Faibussowitsch   PetscCheck(!a->left, PetscObjectComm((PetscObject)C), PETSC_ERR_SUP, "Not implemented");
420fa80d070SPierre Jolivet   B = C->product->B;
4219566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(C, &m, &n));
4229566063dSJacob Faibussowitsch   PetscCall(MatGetSize(C, &M, &N));
423fa80d070SPierre Jolivet   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
4249566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(B, NULL, &n));
4259566063dSJacob Faibussowitsch     PetscCall(MatGetSize(B, NULL, &N));
4269566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &m, NULL));
4279566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &M, NULL));
4289566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, m, n, M, N));
429fa80d070SPierre Jolivet   }
4309566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
4319566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
4329566063dSJacob Faibussowitsch   PetscCall(PetscNew(&contents));
433fa80d070SPierre Jolivet   C->product->data    = contents;
434fa80d070SPierre Jolivet   C->product->destroy = MatNormal_DenseDestroy;
435fa80d070SPierre Jolivet   if (a->right) {
4369566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(a->A, C, NULL, contents->work));
437fa80d070SPierre Jolivet   } else {
4389566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(a->A, B, NULL, contents->work));
439fa80d070SPierre Jolivet   }
4409566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(contents->work[0], MATPRODUCT_AB));
4419566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(contents->work[0]));
4429566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(contents->work[0]));
4439566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(a->A, contents->work[0], NULL, contents->work + 1));
4449566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(contents->work[1], MATPRODUCT_AtB));
4459566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(contents->work[1]));
4469566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(contents->work[1]));
4479566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
4489566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(contents->work[1], array));
4499566063dSJacob Faibussowitsch   PetscCall(MatMPIDenseSetPreallocation(contents->work[1], array));
4509566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
451fa80d070SPierre Jolivet   C->ops->productnumeric = MatProductNumeric_Normal_Dense;
4523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
453fa80d070SPierre Jolivet }
454fa80d070SPierre Jolivet 
455d71ae5a4SJacob Faibussowitsch PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C)
456d71ae5a4SJacob Faibussowitsch {
457fa80d070SPierre Jolivet   PetscFunctionBegin;
458fa80d070SPierre Jolivet   C->ops->productsymbolic = MatProductSymbolic_Normal_Dense;
4593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
460fa80d070SPierre Jolivet }
461fa80d070SPierre Jolivet 
462d71ae5a4SJacob Faibussowitsch PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C)
463d71ae5a4SJacob Faibussowitsch {
464fa80d070SPierre Jolivet   Mat_Product *product = C->product;
465fa80d070SPierre Jolivet 
466fa80d070SPierre Jolivet   PetscFunctionBegin;
46748a46eb9SPierre Jolivet   if (product->type == MATPRODUCT_AB) PetscCall(MatProductSetFromOptions_Normal_Dense_AB(C));
4683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
469fa80d070SPierre Jolivet }
470fa80d070SPierre Jolivet 
4712ef1f0ffSBarry Smith /*MC
4722ef1f0ffSBarry Smith   MATNORMAL - a matrix that behaves like A'*A for `MatMult()` while only containing A
4732ef1f0ffSBarry Smith 
4742ef1f0ffSBarry Smith   Level: intermediate
4752ef1f0ffSBarry Smith 
476*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
4772ef1f0ffSBarry Smith M*/
4782ef1f0ffSBarry Smith 
479c8a8475eSBarry Smith /*@
48011a5261eSBarry Smith       MatCreateNormal - Creates a new `MATNORMAL` matrix object that behaves like A'*A.
481c8a8475eSBarry Smith 
482c3339decSBarry Smith    Collective
483c8a8475eSBarry Smith 
484c8a8475eSBarry Smith    Input Parameter:
485c8a8475eSBarry Smith .   A  - the (possibly rectangular) matrix
486c8a8475eSBarry Smith 
487c8a8475eSBarry Smith    Output Parameter:
488c8a8475eSBarry Smith .   N - the matrix that represents A'*A
489c8a8475eSBarry Smith 
490c30e7cdbSSatish Balay    Level: intermediate
491c30e7cdbSSatish Balay 
49295452b02SPatrick Sanan    Notes:
49395452b02SPatrick Sanan     The product A'*A is NOT actually formed! Rather the new matrix
49411a5261eSBarry Smith           object performs the matrix-vector product, `MatMult()`, by first multiplying by
495c8a8475eSBarry Smith           A and then A'
49611a5261eSBarry Smith 
497*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
498c8a8475eSBarry Smith @*/
499d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormal(Mat A, Mat *N)
500d71ae5a4SJacob Faibussowitsch {
5019ca0e1b6SStefano Zampini   PetscInt    n, nn;
502c8a8475eSBarry Smith   Mat_Normal *Na;
5039ca0e1b6SStefano Zampini   VecType     vtype;
504c8a8475eSBarry Smith 
505c8a8475eSBarry Smith   PetscFunctionBegin;
5069566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &nn));
5079566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, NULL, &n));
5089566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N));
5099566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*N, n, n, nn, nn));
5109566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMAL));
5119566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap));
5129566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap));
513c8a8475eSBarry Smith 
5144dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&Na));
51538f2d2fdSLisandro Dalcin   (*N)->data = (void *)Na;
5169566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
517c3122656SLisandro Dalcin   Na->A     = A;
518b20f02adSBarry Smith   Na->scale = 1.0;
519c8a8475eSBarry Smith 
5209566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, NULL, &Na->w));
5212205254eSKarl Rupp 
522c8a8475eSBarry Smith   (*N)->ops->destroy           = MatDestroy_Normal;
523c8a8475eSBarry Smith   (*N)->ops->mult              = MatMult_Normal;
524b20f02adSBarry Smith   (*N)->ops->multtranspose     = MatMultTranspose_Normal;
525b20f02adSBarry Smith   (*N)->ops->multtransposeadd  = MatMultTransposeAdd_Normal;
526db090513SMatthew Knepley   (*N)->ops->multadd           = MatMultAdd_Normal;
52717a6fd94SBarry Smith   (*N)->ops->getdiagonal       = MatGetDiagonal_Normal;
528a48507d3SJose E. Roman   (*N)->ops->getdiagonalblock  = MatGetDiagonalBlock_Normal;
52969ef1043SBarry Smith   (*N)->ops->scale             = MatScale_Normal;
53069ef1043SBarry Smith   (*N)->ops->diagonalscale     = MatDiagonalScale_Normal;
531fa80d070SPierre Jolivet   (*N)->ops->increaseoverlap   = MatIncreaseOverlap_Normal;
532fa80d070SPierre Jolivet   (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal;
533fa80d070SPierre Jolivet   (*N)->ops->permute           = MatPermute_Normal;
534fa80d070SPierre Jolivet   (*N)->ops->duplicate         = MatDuplicate_Normal;
535fa80d070SPierre Jolivet   (*N)->ops->copy              = MatCopy_Normal;
536c8a8475eSBarry Smith   (*N)->assembled              = PETSC_TRUE;
53748331cefSBarry Smith   (*N)->preallocated           = PETSC_TRUE;
5389ca0e1b6SStefano Zampini 
5399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatNormalGetMat_C", MatNormalGetMat_Normal));
5409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normal_seqaij_C", MatConvert_Normal_AIJ));
5419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normal_mpiaij_C", MatConvert_Normal_AIJ));
5422f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
5432f0a5c5aSJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normal_hypre_C", MatConvert_Normal_HYPRE));
5442f0a5c5aSJose E. Roman #endif
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatProductSetFromOptions_normal_seqdense_C", MatProductSetFromOptions_Normal_Dense));
5469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatProductSetFromOptions_normal_mpidense_C", MatProductSetFromOptions_Normal_Dense));
5479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatProductSetFromOptions_normal_dense_C", MatProductSetFromOptions_Normal_Dense));
5489566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*N, MAT_SYMMETRIC, PETSC_TRUE));
5499566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A, &vtype));
5509566063dSJacob Faibussowitsch   PetscCall(MatSetVecType(*N, vtype));
5519ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE)
5529566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(*N, A->boundtocpu));
5539ca0e1b6SStefano Zampini #endif
5543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
555c8a8475eSBarry Smith }
556