xref: /petsc/src/mat/impls/normal/normm.c (revision 8f83a4d96e6432e9a06bd6cd8c5929cf0bb638f8)
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 
119371c9d4SSatish Balay PetscErrorCode MatScale_Normal(Mat inA, PetscScalar scale) {
1269ef1043SBarry Smith   Mat_Normal *a = (Mat_Normal *)inA->data;
1369ef1043SBarry Smith 
1469ef1043SBarry Smith   PetscFunctionBegin;
1569ef1043SBarry Smith   a->scale *= scale;
1669ef1043SBarry Smith   PetscFunctionReturn(0);
1769ef1043SBarry Smith }
1869ef1043SBarry Smith 
199371c9d4SSatish Balay PetscErrorCode MatDiagonalScale_Normal(Mat inA, Vec left, Vec right) {
207cfd0b8cSBarry Smith   Mat_Normal *a = (Mat_Normal *)inA->data;
217cfd0b8cSBarry Smith 
227cfd0b8cSBarry Smith   PetscFunctionBegin;
237cfd0b8cSBarry Smith   if (left) {
247cfd0b8cSBarry Smith     if (!a->left) {
259566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(left, &a->left));
269566063dSJacob Faibussowitsch       PetscCall(VecCopy(left, a->left));
277cfd0b8cSBarry Smith     } else {
289566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->left, left, a->left));
297cfd0b8cSBarry Smith     }
307cfd0b8cSBarry Smith   }
317cfd0b8cSBarry Smith   if (right) {
327cfd0b8cSBarry Smith     if (!a->right) {
339566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(right, &a->right));
349566063dSJacob Faibussowitsch       PetscCall(VecCopy(right, a->right));
357cfd0b8cSBarry Smith     } else {
369566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->right, right, a->right));
377cfd0b8cSBarry Smith     }
387cfd0b8cSBarry Smith   }
397cfd0b8cSBarry Smith   PetscFunctionReturn(0);
407cfd0b8cSBarry Smith }
417cfd0b8cSBarry Smith 
429371c9d4SSatish Balay PetscErrorCode MatIncreaseOverlap_Normal(Mat A, PetscInt is_max, IS is[], PetscInt ov) {
43fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
44fa80d070SPierre Jolivet   Mat         pattern;
45fa80d070SPierre Jolivet 
46fa80d070SPierre Jolivet   PetscFunctionBegin;
4708401ef6SPierre Jolivet   PetscCheck(ov >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified");
489566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(a->A, a->A, NULL, &pattern));
499566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(pattern, MATPRODUCT_AtB));
509566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(pattern));
519566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(pattern));
529566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(pattern, is_max, is, ov));
539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pattern));
54fa80d070SPierre Jolivet   PetscFunctionReturn(0);
55fa80d070SPierre Jolivet }
56fa80d070SPierre Jolivet 
579371c9d4SSatish Balay PetscErrorCode MatCreateSubMatrices_Normal(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) {
58fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)mat->data;
59fa80d070SPierre Jolivet   Mat         B = a->A, *suba;
60fa80d070SPierre Jolivet   IS         *row;
61fa80d070SPierre Jolivet   PetscInt    M;
62fa80d070SPierre Jolivet 
63fa80d070SPierre Jolivet   PetscFunctionBegin;
64aed4548fSBarry Smith   PetscCheck(!a->left && !a->right && irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented");
6548a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
669566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &M, NULL));
679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &row));
689566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0]));
699566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(row[0]));
70fa80d070SPierre Jolivet   for (M = 1; M < n; ++M) row[M] = row[0];
719566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba));
72fa80d070SPierre Jolivet   for (M = 0; M < n; ++M) {
739566063dSJacob Faibussowitsch     PetscCall(MatCreateNormal(suba[M], *submat + M));
74fa80d070SPierre Jolivet     ((Mat_Normal *)(*submat)[M]->data)->scale = a->scale;
75fa80d070SPierre Jolivet   }
769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&row[0]));
779566063dSJacob Faibussowitsch   PetscCall(PetscFree(row));
789566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(n, &suba));
79fa80d070SPierre Jolivet   PetscFunctionReturn(0);
80fa80d070SPierre Jolivet }
81fa80d070SPierre Jolivet 
829371c9d4SSatish Balay PetscErrorCode MatPermute_Normal(Mat A, IS rowp, IS colp, Mat *B) {
83fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
84fa80d070SPierre Jolivet   Mat         C, Aa = a->A;
85fa80d070SPierre Jolivet   IS          row;
86fa80d070SPierre Jolivet 
87fa80d070SPierre Jolivet   PetscFunctionBegin;
8808401ef6SPierre Jolivet   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
899566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row));
909566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(row));
919566063dSJacob Faibussowitsch   PetscCall(MatPermute(Aa, row, colp, &C));
929566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&row));
939566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(C, B));
949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
95fa80d070SPierre Jolivet   PetscFunctionReturn(0);
96fa80d070SPierre Jolivet }
97fa80d070SPierre Jolivet 
989371c9d4SSatish Balay PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B) {
99fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
100fa80d070SPierre Jolivet   Mat         C;
101fa80d070SPierre Jolivet 
102fa80d070SPierre Jolivet   PetscFunctionBegin;
10308401ef6SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
1049566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(a->A, op, &C));
1059566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(C, B));
1069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
107fa80d070SPierre Jolivet   if (op == MAT_COPY_VALUES) ((Mat_Normal *)(*B)->data)->scale = a->scale;
108fa80d070SPierre Jolivet   PetscFunctionReturn(0);
109fa80d070SPierre Jolivet }
110fa80d070SPierre Jolivet 
1119371c9d4SSatish Balay PetscErrorCode MatCopy_Normal(Mat A, Mat B, MatStructure str) {
112fa80d070SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data, *b = (Mat_Normal *)B->data;
113fa80d070SPierre Jolivet 
114fa80d070SPierre Jolivet   PetscFunctionBegin;
11508401ef6SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
1169566063dSJacob Faibussowitsch   PetscCall(MatCopy(a->A, b->A, str));
117fa80d070SPierre Jolivet   b->scale = a->scale;
1189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->left));
1199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->right));
1209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->leftwork));
1219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->rightwork));
122fa80d070SPierre Jolivet   PetscFunctionReturn(0);
123fa80d070SPierre Jolivet }
124fa80d070SPierre Jolivet 
1259371c9d4SSatish Balay PetscErrorCode MatMult_Normal(Mat N, Vec x, Vec y) {
126c8a8475eSBarry Smith   Mat_Normal *Na = (Mat_Normal *)N->data;
1277cfd0b8cSBarry Smith   Vec         in;
128c8a8475eSBarry Smith 
129c8a8475eSBarry Smith   PetscFunctionBegin;
1307cfd0b8cSBarry Smith   in = x;
1317cfd0b8cSBarry Smith   if (Na->right) {
13248a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1339566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
1347cfd0b8cSBarry Smith     in = Na->rightwork;
1357cfd0b8cSBarry Smith   }
1369566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1379566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(Na->A, Na->w, y));
1381baa6e33SBarry Smith   if (Na->left) PetscCall(VecPointwiseMult(y, Na->left, y));
1399566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
140c8a8475eSBarry Smith   PetscFunctionReturn(0);
141c8a8475eSBarry Smith }
142c8a8475eSBarry Smith 
1439371c9d4SSatish Balay PetscErrorCode MatMultAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3) {
144db090513SMatthew Knepley   Mat_Normal *Na = (Mat_Normal *)N->data;
1457cfd0b8cSBarry Smith   Vec         in;
146db090513SMatthew Knepley 
147db090513SMatthew Knepley   PetscFunctionBegin;
1487cfd0b8cSBarry Smith   in = v1;
1497cfd0b8cSBarry Smith   if (Na->right) {
15048a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1519566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
1527cfd0b8cSBarry Smith     in = Na->rightwork;
1537cfd0b8cSBarry Smith   }
1549566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1559566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
1567cfd0b8cSBarry Smith   if (Na->left) {
1579566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(Na->A, Na->w, v3));
1589566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->left, v3));
1599566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
1607cfd0b8cSBarry Smith   } else {
1619566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(Na->A, Na->w, v2, v3));
1627cfd0b8cSBarry Smith   }
163b20f02adSBarry Smith   PetscFunctionReturn(0);
164b20f02adSBarry Smith }
165b20f02adSBarry Smith 
1669371c9d4SSatish Balay PetscErrorCode MatMultTranspose_Normal(Mat N, Vec x, Vec y) {
167b20f02adSBarry Smith   Mat_Normal *Na = (Mat_Normal *)N->data;
1687cfd0b8cSBarry Smith   Vec         in;
169b20f02adSBarry Smith 
170b20f02adSBarry Smith   PetscFunctionBegin;
1717cfd0b8cSBarry Smith   in = x;
1727cfd0b8cSBarry Smith   if (Na->left) {
17348a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
1749566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
1757cfd0b8cSBarry Smith     in = Na->leftwork;
1767cfd0b8cSBarry Smith   }
1779566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1789566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(Na->A, Na->w, y));
1791baa6e33SBarry Smith   if (Na->right) PetscCall(VecPointwiseMult(y, Na->right, y));
1809566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
181b20f02adSBarry Smith   PetscFunctionReturn(0);
182b20f02adSBarry Smith }
183b20f02adSBarry Smith 
1849371c9d4SSatish Balay PetscErrorCode MatMultTransposeAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3) {
185b20f02adSBarry Smith   Mat_Normal *Na = (Mat_Normal *)N->data;
1867cfd0b8cSBarry Smith   Vec         in;
187b20f02adSBarry Smith 
188b20f02adSBarry Smith   PetscFunctionBegin;
1897cfd0b8cSBarry Smith   in = v1;
1907cfd0b8cSBarry Smith   if (Na->left) {
19148a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
1929566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
1937cfd0b8cSBarry Smith     in = Na->leftwork;
1947cfd0b8cSBarry Smith   }
1959566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1969566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
1977cfd0b8cSBarry Smith   if (Na->right) {
1989566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(Na->A, Na->w, v3));
1999566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->right, v3));
2009566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
2017cfd0b8cSBarry Smith   } else {
2029566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(Na->A, Na->w, v2, v3));
2037cfd0b8cSBarry Smith   }
204db090513SMatthew Knepley   PetscFunctionReturn(0);
205db090513SMatthew Knepley }
206db090513SMatthew Knepley 
2079371c9d4SSatish Balay PetscErrorCode MatDestroy_Normal(Mat N) {
208c8a8475eSBarry Smith   Mat_Normal *Na = (Mat_Normal *)N->data;
209c8a8475eSBarry Smith 
210c8a8475eSBarry Smith   PetscFunctionBegin;
2119566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
212a48507d3SJose E. Roman   PetscCall(MatDestroy(&Na->D));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->w));
2149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->left));
2159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->right));
2169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->leftwork));
2179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->rightwork));
2189566063dSJacob Faibussowitsch   PetscCall(PetscFree(N->data));
2199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMat_C", NULL));
2209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_seqaij_C", NULL));
2219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_mpiaij_C", NULL));
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_seqdense_C", NULL));
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_mpidense_C", NULL));
2249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_dense_C", NULL));
225c8a8475eSBarry Smith   PetscFunctionReturn(0);
226c8a8475eSBarry Smith }
227c8a8475eSBarry Smith 
22817a6fd94SBarry Smith /*
22917a6fd94SBarry Smith       Slow, nonscalable version
23017a6fd94SBarry Smith */
2319371c9d4SSatish Balay PetscErrorCode MatGetDiagonal_Normal(Mat N, Vec v) {
23217a6fd94SBarry Smith   Mat_Normal        *Na = (Mat_Normal *)N->data;
23317a6fd94SBarry Smith   Mat                A  = Na->A;
234b24ad042SBarry Smith   PetscInt           i, j, rstart, rend, nnz;
235b24ad042SBarry Smith   const PetscInt    *cols;
23617a6fd94SBarry Smith   PetscScalar       *diag, *work, *values;
237b3cc6726SBarry Smith   const PetscScalar *mvalues;
23817a6fd94SBarry Smith 
23917a6fd94SBarry Smith   PetscFunctionBegin;
2409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work));
2419566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(work, A->cmap->N));
2429566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
24317a6fd94SBarry Smith   for (i = rstart; i < rend; i++) {
2449566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues));
245ad540459SPierre Jolivet     for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * mvalues[j];
2469566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues));
24717a6fd94SBarry Smith   }
2481c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N)));
249d0f46423SBarry Smith   rstart = N->cmap->rstart;
250d0f46423SBarry Smith   rend   = N->cmap->rend;
2519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &values));
2529566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart));
2539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &values));
2549566063dSJacob Faibussowitsch   PetscCall(PetscFree2(diag, work));
2559566063dSJacob Faibussowitsch   PetscCall(VecScale(v, Na->scale));
25617a6fd94SBarry Smith   PetscFunctionReturn(0);
25717a6fd94SBarry Smith }
258c8a8475eSBarry Smith 
2599371c9d4SSatish Balay PetscErrorCode MatGetDiagonalBlock_Normal(Mat N, Mat *D) {
260a48507d3SJose E. Roman   Mat_Normal *Na = (Mat_Normal *)N->data;
261a48507d3SJose E. Roman   Mat         M, A = Na->A;
262a48507d3SJose E. Roman 
263a48507d3SJose E. Roman   PetscFunctionBegin;
264a48507d3SJose E. Roman   PetscCall(MatGetDiagonalBlock(A, &M));
265a48507d3SJose E. Roman   PetscCall(MatCreateNormal(M, &Na->D));
266a48507d3SJose E. Roman   *D = Na->D;
267a48507d3SJose E. Roman   PetscFunctionReturn(0);
268a48507d3SJose E. Roman }
269a48507d3SJose E. Roman 
2709371c9d4SSatish Balay PetscErrorCode MatNormalGetMat_Normal(Mat A, Mat *M) {
271fa80d070SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal *)A->data;
272fa80d070SPierre Jolivet 
273fa80d070SPierre Jolivet   PetscFunctionBegin;
274fa80d070SPierre Jolivet   *M = Aa->A;
275fa80d070SPierre Jolivet   PetscFunctionReturn(0);
276fa80d070SPierre Jolivet }
277fa80d070SPierre Jolivet 
278fa80d070SPierre Jolivet /*@
279fa80d070SPierre Jolivet       MatNormalGetMat - Gets the Mat object stored inside a MATNORMAL
280fa80d070SPierre Jolivet 
281fa80d070SPierre Jolivet    Logically collective on Mat
282fa80d070SPierre Jolivet 
283fa80d070SPierre Jolivet    Input Parameter:
284fa80d070SPierre Jolivet .   A  - the MATNORMAL matrix
285fa80d070SPierre Jolivet 
286fa80d070SPierre Jolivet    Output Parameter:
287fa80d070SPierre Jolivet .   M - the matrix object stored inside A
288fa80d070SPierre Jolivet 
289fa80d070SPierre Jolivet    Level: intermediate
290fa80d070SPierre Jolivet 
291db781477SPatrick Sanan .seealso: `MatCreateNormal()`
292fa80d070SPierre Jolivet 
293fa80d070SPierre Jolivet @*/
2949371c9d4SSatish Balay PetscErrorCode MatNormalGetMat(Mat A, Mat *M) {
295fa80d070SPierre Jolivet   PetscFunctionBegin;
296fa80d070SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
297fa80d070SPierre Jolivet   PetscValidType(A, 1);
298fa80d070SPierre Jolivet   PetscValidPointer(M, 2);
299cac4c232SBarry Smith   PetscUseMethod(A, "MatNormalGetMat_C", (Mat, Mat *), (A, M));
300fa80d070SPierre Jolivet   PetscFunctionReturn(0);
301fa80d070SPierre Jolivet }
302fa80d070SPierre Jolivet 
3039371c9d4SSatish Balay PetscErrorCode MatConvert_Normal_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) {
304fa80d070SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal *)A->data;
305fa80d070SPierre Jolivet   Mat         B;
306fa80d070SPierre Jolivet   PetscInt    m, n, M, N;
307fa80d070SPierre Jolivet 
308fa80d070SPierre Jolivet   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
3109566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
311fa80d070SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
312fa80d070SPierre Jolivet     B = *newmat;
3139566063dSJacob Faibussowitsch     PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B));
314fa80d070SPierre Jolivet   } else {
3159566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B));
3169566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
3179566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(B));
3189566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(B));
3199566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE));
320fa80d070SPierre Jolivet   }
3219566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(B));
322fa80d070SPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
3239566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
324fa80d070SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
3259566063dSJacob Faibussowitsch   PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat));
326fa80d070SPierre Jolivet   PetscFunctionReturn(0);
327fa80d070SPierre Jolivet }
328fa80d070SPierre Jolivet 
329fa80d070SPierre Jolivet typedef struct {
330fa80d070SPierre Jolivet   Mat work[2];
331fa80d070SPierre Jolivet } Normal_Dense;
332fa80d070SPierre Jolivet 
3339371c9d4SSatish Balay PetscErrorCode MatProductNumeric_Normal_Dense(Mat C) {
334fa80d070SPierre Jolivet   Mat           A, B;
335fa80d070SPierre Jolivet   Normal_Dense *contents;
336fa80d070SPierre Jolivet   Mat_Normal   *a;
337fa80d070SPierre Jolivet   PetscScalar  *array;
338fa80d070SPierre Jolivet 
339fa80d070SPierre Jolivet   PetscFunctionBegin;
340*8f83a4d9SJacob Faibussowitsch   MatCheckProduct(C, 1);
341fa80d070SPierre Jolivet   A        = C->product->A;
342fa80d070SPierre Jolivet   a        = (Mat_Normal *)A->data;
343fa80d070SPierre Jolivet   B        = C->product->B;
344fa80d070SPierre Jolivet   contents = (Normal_Dense *)C->product->data;
34528b400f6SJacob Faibussowitsch   PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
346fa80d070SPierre Jolivet   if (a->right) {
3479566063dSJacob Faibussowitsch     PetscCall(MatCopy(B, C, SAME_NONZERO_PATTERN));
3489566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(C, a->right, NULL));
349fa80d070SPierre Jolivet   }
3509566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(contents->work[0]));
3519566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
3529566063dSJacob Faibussowitsch   PetscCall(MatDensePlaceArray(contents->work[1], array));
3539566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(contents->work[1]));
3549566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
3559566063dSJacob Faibussowitsch   PetscCall(MatDenseResetArray(contents->work[1]));
3569566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
3579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
3589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
3599566063dSJacob Faibussowitsch   PetscCall(MatScale(C, a->scale));
360fa80d070SPierre Jolivet   PetscFunctionReturn(0);
361fa80d070SPierre Jolivet }
362fa80d070SPierre Jolivet 
3639371c9d4SSatish Balay PetscErrorCode MatNormal_DenseDestroy(void *ctx) {
364fa80d070SPierre Jolivet   Normal_Dense *contents = (Normal_Dense *)ctx;
365fa80d070SPierre Jolivet 
366fa80d070SPierre Jolivet   PetscFunctionBegin;
3679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(contents->work));
3689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(contents->work + 1));
3699566063dSJacob Faibussowitsch   PetscCall(PetscFree(contents));
370fa80d070SPierre Jolivet   PetscFunctionReturn(0);
371fa80d070SPierre Jolivet }
372fa80d070SPierre Jolivet 
3739371c9d4SSatish Balay PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C) {
374fa80d070SPierre Jolivet   Mat           A, B;
375fa80d070SPierre Jolivet   Normal_Dense *contents = NULL;
376fa80d070SPierre Jolivet   Mat_Normal   *a;
377fa80d070SPierre Jolivet   PetscScalar  *array;
378fa80d070SPierre Jolivet   PetscInt      n, N, m, M;
379fa80d070SPierre Jolivet 
380fa80d070SPierre Jolivet   PetscFunctionBegin;
381*8f83a4d9SJacob Faibussowitsch   MatCheckProduct(C, 1);
38228b400f6SJacob Faibussowitsch   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
383fa80d070SPierre Jolivet   A = C->product->A;
384fa80d070SPierre Jolivet   a = (Mat_Normal *)A->data;
38528b400f6SJacob Faibussowitsch   PetscCheck(!a->left, PetscObjectComm((PetscObject)C), PETSC_ERR_SUP, "Not implemented");
386fa80d070SPierre Jolivet   B = C->product->B;
3879566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(C, &m, &n));
3889566063dSJacob Faibussowitsch   PetscCall(MatGetSize(C, &M, &N));
389fa80d070SPierre Jolivet   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
3909566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(B, NULL, &n));
3919566063dSJacob Faibussowitsch     PetscCall(MatGetSize(B, NULL, &N));
3929566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &m, NULL));
3939566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &M, NULL));
3949566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, m, n, M, N));
395fa80d070SPierre Jolivet   }
3969566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
3979566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
3989566063dSJacob Faibussowitsch   PetscCall(PetscNew(&contents));
399fa80d070SPierre Jolivet   C->product->data    = contents;
400fa80d070SPierre Jolivet   C->product->destroy = MatNormal_DenseDestroy;
401fa80d070SPierre Jolivet   if (a->right) {
4029566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(a->A, C, NULL, contents->work));
403fa80d070SPierre Jolivet   } else {
4049566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(a->A, B, NULL, contents->work));
405fa80d070SPierre Jolivet   }
4069566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(contents->work[0], MATPRODUCT_AB));
4079566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(contents->work[0]));
4089566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(contents->work[0]));
4099566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(a->A, contents->work[0], NULL, contents->work + 1));
4109566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(contents->work[1], MATPRODUCT_AtB));
4119566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(contents->work[1]));
4129566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(contents->work[1]));
4139566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
4149566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(contents->work[1], array));
4159566063dSJacob Faibussowitsch   PetscCall(MatMPIDenseSetPreallocation(contents->work[1], array));
4169566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
417fa80d070SPierre Jolivet   C->ops->productnumeric = MatProductNumeric_Normal_Dense;
418fa80d070SPierre Jolivet   PetscFunctionReturn(0);
419fa80d070SPierre Jolivet }
420fa80d070SPierre Jolivet 
4219371c9d4SSatish Balay PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C) {
422fa80d070SPierre Jolivet   PetscFunctionBegin;
423fa80d070SPierre Jolivet   C->ops->productsymbolic = MatProductSymbolic_Normal_Dense;
424fa80d070SPierre Jolivet   PetscFunctionReturn(0);
425fa80d070SPierre Jolivet }
426fa80d070SPierre Jolivet 
4279371c9d4SSatish Balay PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C) {
428fa80d070SPierre Jolivet   Mat_Product *product = C->product;
429fa80d070SPierre Jolivet 
430fa80d070SPierre Jolivet   PetscFunctionBegin;
43148a46eb9SPierre Jolivet   if (product->type == MATPRODUCT_AB) PetscCall(MatProductSetFromOptions_Normal_Dense_AB(C));
432fa80d070SPierre Jolivet   PetscFunctionReturn(0);
433fa80d070SPierre Jolivet }
434fa80d070SPierre Jolivet 
435c8a8475eSBarry Smith /*@
436c8a8475eSBarry Smith       MatCreateNormal - Creates a new matrix object that behaves like A'*A.
437c8a8475eSBarry Smith 
438c8a8475eSBarry Smith    Collective on Mat
439c8a8475eSBarry Smith 
440c8a8475eSBarry Smith    Input Parameter:
441c8a8475eSBarry Smith .   A  - the (possibly rectangular) matrix
442c8a8475eSBarry Smith 
443c8a8475eSBarry Smith    Output Parameter:
444c8a8475eSBarry Smith .   N - the matrix that represents A'*A
445c8a8475eSBarry Smith 
446c30e7cdbSSatish Balay    Level: intermediate
447c30e7cdbSSatish Balay 
44895452b02SPatrick Sanan    Notes:
44995452b02SPatrick Sanan     The product A'*A is NOT actually formed! Rather the new matrix
450c8a8475eSBarry Smith           object performs the matrix-vector product by first multiplying by
451c8a8475eSBarry Smith           A and then A'
452c8a8475eSBarry Smith @*/
4539371c9d4SSatish Balay PetscErrorCode MatCreateNormal(Mat A, Mat *N) {
4549ca0e1b6SStefano Zampini   PetscInt    n, nn;
455c8a8475eSBarry Smith   Mat_Normal *Na;
4569ca0e1b6SStefano Zampini   VecType     vtype;
457c8a8475eSBarry Smith 
458c8a8475eSBarry Smith   PetscFunctionBegin;
4599566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &nn));
4609566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, NULL, &n));
4619566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N));
4629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*N, n, n, nn, nn));
4639566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMAL));
4649566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap));
4659566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap));
466c8a8475eSBarry Smith 
4679566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(*N, &Na));
46838f2d2fdSLisandro Dalcin   (*N)->data = (void *)Na;
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
470c3122656SLisandro Dalcin   Na->A     = A;
471b20f02adSBarry Smith   Na->scale = 1.0;
472c8a8475eSBarry Smith 
4739566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, NULL, &Na->w));
4742205254eSKarl Rupp 
475c8a8475eSBarry Smith   (*N)->ops->destroy           = MatDestroy_Normal;
476c8a8475eSBarry Smith   (*N)->ops->mult              = MatMult_Normal;
477b20f02adSBarry Smith   (*N)->ops->multtranspose     = MatMultTranspose_Normal;
478b20f02adSBarry Smith   (*N)->ops->multtransposeadd  = MatMultTransposeAdd_Normal;
479db090513SMatthew Knepley   (*N)->ops->multadd           = MatMultAdd_Normal;
48017a6fd94SBarry Smith   (*N)->ops->getdiagonal       = MatGetDiagonal_Normal;
481a48507d3SJose E. Roman   (*N)->ops->getdiagonalblock  = MatGetDiagonalBlock_Normal;
48269ef1043SBarry Smith   (*N)->ops->scale             = MatScale_Normal;
48369ef1043SBarry Smith   (*N)->ops->diagonalscale     = MatDiagonalScale_Normal;
484fa80d070SPierre Jolivet   (*N)->ops->increaseoverlap   = MatIncreaseOverlap_Normal;
485fa80d070SPierre Jolivet   (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal;
486fa80d070SPierre Jolivet   (*N)->ops->permute           = MatPermute_Normal;
487fa80d070SPierre Jolivet   (*N)->ops->duplicate         = MatDuplicate_Normal;
488fa80d070SPierre Jolivet   (*N)->ops->copy              = MatCopy_Normal;
489c8a8475eSBarry Smith   (*N)->assembled              = PETSC_TRUE;
49048331cefSBarry Smith   (*N)->preallocated           = PETSC_TRUE;
4919ca0e1b6SStefano Zampini 
4929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatNormalGetMat_C", MatNormalGetMat_Normal));
4939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normal_seqaij_C", MatConvert_Normal_AIJ));
4949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normal_mpiaij_C", MatConvert_Normal_AIJ));
4959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatProductSetFromOptions_normal_seqdense_C", MatProductSetFromOptions_Normal_Dense));
4969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatProductSetFromOptions_normal_mpidense_C", MatProductSetFromOptions_Normal_Dense));
4979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatProductSetFromOptions_normal_dense_C", MatProductSetFromOptions_Normal_Dense));
4989566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*N, MAT_SYMMETRIC, PETSC_TRUE));
4999566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A, &vtype));
5009566063dSJacob Faibussowitsch   PetscCall(MatSetVecType(*N, vtype));
5019ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE)
5029566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(*N, A->boundtocpu));
5039ca0e1b6SStefano Zampini #endif
504c8a8475eSBarry Smith   PetscFunctionReturn(0);
505c8a8475eSBarry Smith }
506