xref: /petsc/src/mat/impls/normal/normmh.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
1ebda224cSfranklin5 
2ebda224cSfranklin5 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
3ebda224cSfranklin5 
4ebda224cSfranklin5 typedef struct {
5ebda224cSfranklin5   Mat         A;
6a48507d3SJose E. Roman   Mat         D; /* local submatrix for diagonal part */
7ebda224cSfranklin5   Vec         w, left, right, leftwork, rightwork;
8ebda224cSfranklin5   PetscScalar scale;
9ebda224cSfranklin5 } Mat_Normal;
10ebda224cSfranklin5 
11d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_NormalHermitian(Mat inA, PetscScalar scale)
12d71ae5a4SJacob Faibussowitsch {
13ebda224cSfranklin5   Mat_Normal *a = (Mat_Normal *)inA->data;
14ebda224cSfranklin5 
15ebda224cSfranklin5   PetscFunctionBegin;
16ebda224cSfranklin5   a->scale *= scale;
173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18ebda224cSfranklin5 }
19ebda224cSfranklin5 
20d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalScale_NormalHermitian(Mat inA, Vec left, Vec right)
21d71ae5a4SJacob Faibussowitsch {
22ebda224cSfranklin5   Mat_Normal *a = (Mat_Normal *)inA->data;
23ebda224cSfranklin5 
24ebda224cSfranklin5   PetscFunctionBegin;
25ebda224cSfranklin5   if (left) {
26ebda224cSfranklin5     if (!a->left) {
279566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(left, &a->left));
289566063dSJacob Faibussowitsch       PetscCall(VecCopy(left, a->left));
29ebda224cSfranklin5     } else {
309566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->left, left, a->left));
31ebda224cSfranklin5     }
32ebda224cSfranklin5   }
33ebda224cSfranklin5   if (right) {
34ebda224cSfranklin5     if (!a->right) {
359566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(right, &a->right));
369566063dSJacob Faibussowitsch       PetscCall(VecCopy(right, a->right));
37ebda224cSfranklin5     } else {
389566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->right, right, a->right));
39ebda224cSfranklin5     }
40ebda224cSfranklin5   }
413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42ebda224cSfranklin5 }
43ebda224cSfranklin5 
44d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrices_NormalHermitian(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
45d71ae5a4SJacob Faibussowitsch {
46ab704866SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)mat->data;
47ab704866SPierre Jolivet   Mat         B = a->A, *suba;
48ab704866SPierre Jolivet   IS         *row;
49ab704866SPierre Jolivet   PetscInt    M;
50ab704866SPierre Jolivet 
51ab704866SPierre Jolivet   PetscFunctionBegin;
52ab704866SPierre Jolivet   PetscCheck(!a->left && !a->right && irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented");
5348a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
54ab704866SPierre Jolivet   PetscCall(MatGetSize(B, &M, NULL));
55ab704866SPierre Jolivet   PetscCall(PetscMalloc1(n, &row));
56ab704866SPierre Jolivet   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0]));
57ab704866SPierre Jolivet   PetscCall(ISSetIdentity(row[0]));
58ab704866SPierre Jolivet   for (M = 1; M < n; ++M) row[M] = row[0];
59ab704866SPierre Jolivet   PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba));
60ab704866SPierre Jolivet   for (M = 0; M < n; ++M) {
61ab704866SPierre Jolivet     PetscCall(MatCreateNormalHermitian(suba[M], *submat + M));
62ab704866SPierre Jolivet     ((Mat_Normal *)(*submat)[M]->data)->scale = a->scale;
63ab704866SPierre Jolivet   }
64ab704866SPierre Jolivet   PetscCall(ISDestroy(&row[0]));
65ab704866SPierre Jolivet   PetscCall(PetscFree(row));
66ab704866SPierre Jolivet   PetscCall(MatDestroySubMatrices(n, &suba));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68ab704866SPierre Jolivet }
69ab704866SPierre Jolivet 
70d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPermute_NormalHermitian(Mat A, IS rowp, IS colp, Mat *B)
71d71ae5a4SJacob Faibussowitsch {
72ab704866SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
73ab704866SPierre Jolivet   Mat         C, Aa = a->A;
74ab704866SPierre Jolivet   IS          row;
75ab704866SPierre Jolivet 
76ab704866SPierre Jolivet   PetscFunctionBegin;
77ab704866SPierre Jolivet   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
78ab704866SPierre Jolivet   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row));
79ab704866SPierre Jolivet   PetscCall(ISSetIdentity(row));
80ab704866SPierre Jolivet   PetscCall(MatPermute(Aa, row, colp, &C));
81ab704866SPierre Jolivet   PetscCall(ISDestroy(&row));
82ab704866SPierre Jolivet   PetscCall(MatCreateNormalHermitian(C, B));
83ab704866SPierre Jolivet   PetscCall(MatDestroy(&C));
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
85ab704866SPierre Jolivet }
86ab704866SPierre Jolivet 
87d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_NormalHermitian(Mat A, MatDuplicateOption op, Mat *B)
88d71ae5a4SJacob Faibussowitsch {
89ab704866SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
90ab704866SPierre Jolivet   Mat         C;
91ab704866SPierre Jolivet 
92ab704866SPierre Jolivet   PetscFunctionBegin;
93ab704866SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
94ab704866SPierre Jolivet   PetscCall(MatDuplicate(a->A, op, &C));
95ab704866SPierre Jolivet   PetscCall(MatCreateNormalHermitian(C, B));
96ab704866SPierre Jolivet   PetscCall(MatDestroy(&C));
97ab704866SPierre Jolivet   if (op == MAT_COPY_VALUES) ((Mat_Normal *)(*B)->data)->scale = a->scale;
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99ab704866SPierre Jolivet }
100ab704866SPierre Jolivet 
101d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCopy_NormalHermitian(Mat A, Mat B, MatStructure str)
102d71ae5a4SJacob Faibussowitsch {
103ab704866SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data, *b = (Mat_Normal *)B->data;
104ab704866SPierre Jolivet 
105ab704866SPierre Jolivet   PetscFunctionBegin;
106ab704866SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
107ab704866SPierre Jolivet   PetscCall(MatCopy(a->A, b->A, str));
108ab704866SPierre Jolivet   b->scale = a->scale;
109ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->left));
110ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->right));
111ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->leftwork));
112ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->rightwork));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114ab704866SPierre Jolivet }
115ab704866SPierre Jolivet 
116d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_NormalHermitian(Mat N, Vec x, Vec y)
117d71ae5a4SJacob Faibussowitsch {
118ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
119ebda224cSfranklin5   Vec         in;
120ebda224cSfranklin5 
121ebda224cSfranklin5   PetscFunctionBegin;
122ebda224cSfranklin5   in = x;
123ebda224cSfranklin5   if (Na->right) {
12448a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1259566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
126ebda224cSfranklin5     in = Na->rightwork;
127ebda224cSfranklin5   }
1289566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1299566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y));
1301baa6e33SBarry Smith   if (Na->left) PetscCall(VecPointwiseMult(y, Na->left, y));
1319566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
1323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
133ebda224cSfranklin5 }
134ebda224cSfranklin5 
135d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3)
136d71ae5a4SJacob Faibussowitsch {
137ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
138ebda224cSfranklin5   Vec         in;
139ebda224cSfranklin5 
140ebda224cSfranklin5   PetscFunctionBegin;
141ebda224cSfranklin5   in = v1;
142ebda224cSfranklin5   if (Na->right) {
14348a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1449566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
145ebda224cSfranklin5     in = Na->rightwork;
146ebda224cSfranklin5   }
1479566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1489566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
149ebda224cSfranklin5   if (Na->left) {
1509566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(Na->A, Na->w, v3));
1519566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->left, v3));
1529566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
153ebda224cSfranklin5   } else {
1549566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTransposeAdd(Na->A, Na->w, v2, v3));
155ebda224cSfranklin5   }
1563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
157ebda224cSfranklin5 }
158ebda224cSfranklin5 
159d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTranspose_Normal(Mat N, Vec x, Vec y)
160d71ae5a4SJacob Faibussowitsch {
161ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
162ebda224cSfranklin5   Vec         in;
163ebda224cSfranklin5 
164ebda224cSfranklin5   PetscFunctionBegin;
165ebda224cSfranklin5   in = x;
166ebda224cSfranklin5   if (Na->left) {
16748a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
1689566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
169ebda224cSfranklin5     in = Na->leftwork;
170ebda224cSfranklin5   }
1719566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1729566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y));
1731baa6e33SBarry Smith   if (Na->right) PetscCall(VecPointwiseMult(y, Na->right, y));
1749566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
1753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
176ebda224cSfranklin5 }
177ebda224cSfranklin5 
178d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3)
179d71ae5a4SJacob Faibussowitsch {
180ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
181ebda224cSfranklin5   Vec         in;
182ebda224cSfranklin5 
183ebda224cSfranklin5   PetscFunctionBegin;
184ebda224cSfranklin5   in = v1;
185ebda224cSfranklin5   if (Na->left) {
18648a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
1879566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
188ebda224cSfranklin5     in = Na->leftwork;
189ebda224cSfranklin5   }
1909566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1919566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
192ebda224cSfranklin5   if (Na->right) {
1939566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(Na->A, Na->w, v3));
1949566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->right, v3));
1959566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
196ebda224cSfranklin5   } else {
1979566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTransposeAdd(Na->A, Na->w, v2, v3));
198ebda224cSfranklin5   }
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
200ebda224cSfranklin5 }
201ebda224cSfranklin5 
202d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_NormalHermitian(Mat N)
203d71ae5a4SJacob Faibussowitsch {
204ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
205ebda224cSfranklin5 
206ebda224cSfranklin5   PetscFunctionBegin;
2079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
208a48507d3SJose E. Roman   PetscCall(MatDestroy(&Na->D));
2099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->w));
2109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->left));
2119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->right));
2129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->leftwork));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->rightwork));
2149566063dSJacob Faibussowitsch   PetscCall(PetscFree(N->data));
2159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMatHermitian_C", NULL));
216ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_seqaij_C", NULL));
217ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_mpiaij_C", NULL));
2182f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
2192f0a5c5aSJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_hypre_C", NULL));
2202f0a5c5aSJose E. Roman #endif
2213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
222ebda224cSfranklin5 }
223ebda224cSfranklin5 
224ebda224cSfranklin5 /*
225ebda224cSfranklin5       Slow, nonscalable version
226ebda224cSfranklin5 */
227d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_NormalHermitian(Mat N, Vec v)
228d71ae5a4SJacob Faibussowitsch {
229ebda224cSfranklin5   Mat_Normal        *Na = (Mat_Normal *)N->data;
230ebda224cSfranklin5   Mat                A  = Na->A;
231ebda224cSfranklin5   PetscInt           i, j, rstart, rend, nnz;
232ebda224cSfranklin5   const PetscInt    *cols;
233ebda224cSfranklin5   PetscScalar       *diag, *work, *values;
234ebda224cSfranklin5   const PetscScalar *mvalues;
235ebda224cSfranklin5 
236ebda224cSfranklin5   PetscFunctionBegin;
2379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work));
2389566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(work, A->cmap->N));
2399566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
240ebda224cSfranklin5   for (i = rstart; i < rend; i++) {
2419566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues));
242ad540459SPierre Jolivet     for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * PetscConj(mvalues[j]);
2439566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues));
244ebda224cSfranklin5   }
2451c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N)));
246ebda224cSfranklin5   rstart = N->cmap->rstart;
247ebda224cSfranklin5   rend   = N->cmap->rend;
2489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &values));
2499566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart));
2509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &values));
2519566063dSJacob Faibussowitsch   PetscCall(PetscFree2(diag, work));
2529566063dSJacob Faibussowitsch   PetscCall(VecScale(v, Na->scale));
2533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
254ebda224cSfranklin5 }
255ebda224cSfranklin5 
256d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonalBlock_NormalHermitian(Mat N, Mat *D)
257d71ae5a4SJacob Faibussowitsch {
258a48507d3SJose E. Roman   Mat_Normal *Na = (Mat_Normal *)N->data;
259a48507d3SJose E. Roman   Mat         M, A = Na->A;
260a48507d3SJose E. Roman 
261a48507d3SJose E. Roman   PetscFunctionBegin;
262a48507d3SJose E. Roman   PetscCall(MatGetDiagonalBlock(A, &M));
263a48507d3SJose E. Roman   PetscCall(MatCreateNormalHermitian(M, &Na->D));
264a48507d3SJose E. Roman   *D = Na->D;
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
266a48507d3SJose E. Roman }
267a48507d3SJose E. Roman 
268d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalGetMat_NormalHermitian(Mat A, Mat *M)
269d71ae5a4SJacob Faibussowitsch {
27065f45395SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal *)A->data;
27165f45395SPierre Jolivet 
27265f45395SPierre Jolivet   PetscFunctionBegin;
27365f45395SPierre Jolivet   *M = Aa->A;
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27565f45395SPierre Jolivet }
27665f45395SPierre Jolivet 
27765f45395SPierre Jolivet /*@
27811a5261eSBarry Smith       MatNormalHermitianGetMat - Gets the `Mat` object stored inside a `MATNORMALHERMITIAN`
27965f45395SPierre Jolivet 
2802ef1f0ffSBarry Smith    Logically Collective
28165f45395SPierre Jolivet 
28265f45395SPierre Jolivet    Input Parameter:
28311a5261eSBarry Smith .   A  - the `MATNORMALHERMITIAN` matrix
28465f45395SPierre Jolivet 
28565f45395SPierre Jolivet    Output Parameter:
28665f45395SPierre Jolivet .   M - the matrix object stored inside A
28765f45395SPierre Jolivet 
28865f45395SPierre Jolivet    Level: intermediate
28965f45395SPierre Jolivet 
290*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
29165f45395SPierre Jolivet @*/
292d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalHermitianGetMat(Mat A, Mat *M)
293d71ae5a4SJacob Faibussowitsch {
29465f45395SPierre Jolivet   PetscFunctionBegin;
29565f45395SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
29665f45395SPierre Jolivet   PetscValidType(A, 1);
29765f45395SPierre Jolivet   PetscValidPointer(M, 2);
298cac4c232SBarry Smith   PetscUseMethod(A, "MatNormalGetMatHermitian_C", (Mat, Mat *), (A, M));
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30065f45395SPierre Jolivet }
30165f45395SPierre Jolivet 
302d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_NormalHermitian_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
303d71ae5a4SJacob Faibussowitsch {
304ab704866SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal *)A->data;
305ab704866SPierre Jolivet   Mat         B, conjugate;
306ab704866SPierre Jolivet   PetscInt    m, n, M, N;
307ab704866SPierre Jolivet 
308ab704866SPierre Jolivet   PetscFunctionBegin;
309ab704866SPierre Jolivet   PetscCall(MatGetSize(A, &M, &N));
310ab704866SPierre Jolivet   PetscCall(MatGetLocalSize(A, &m, &n));
311ab704866SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
312ab704866SPierre Jolivet     B = *newmat;
313ab704866SPierre Jolivet     PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B));
314ab704866SPierre Jolivet   } else {
315ab704866SPierre Jolivet     PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B));
316ab704866SPierre Jolivet     PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
317ab704866SPierre Jolivet     PetscCall(MatProductSetFromOptions(B));
318ab704866SPierre Jolivet     PetscCall(MatProductSymbolic(B));
319ab704866SPierre Jolivet     PetscCall(MatSetOption(B, !PetscDefined(USE_COMPLEX) ? MAT_SYMMETRIC : MAT_HERMITIAN, PETSC_TRUE));
320ab704866SPierre Jolivet   }
321ab704866SPierre Jolivet   if (PetscDefined(USE_COMPLEX)) {
322ab704866SPierre Jolivet     PetscCall(MatDuplicate(Aa->A, MAT_COPY_VALUES, &conjugate));
323ab704866SPierre Jolivet     PetscCall(MatConjugate(conjugate));
324ab704866SPierre Jolivet     PetscCall(MatProductReplaceMats(conjugate, Aa->A, NULL, B));
325ab704866SPierre Jolivet   }
326ab704866SPierre Jolivet   PetscCall(MatProductNumeric(B));
327ab704866SPierre Jolivet   if (PetscDefined(USE_COMPLEX)) PetscCall(MatDestroy(&conjugate));
328ab704866SPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
329ab704866SPierre Jolivet     PetscCall(MatHeaderReplace(A, &B));
330ab704866SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
331ab704866SPierre Jolivet   PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat));
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
333ab704866SPierre Jolivet }
334ab704866SPierre Jolivet 
3352f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
3362f0a5c5aSJose E. Roman PetscErrorCode MatConvert_NormalHermitian_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
3372f0a5c5aSJose E. Roman {
3382f0a5c5aSJose E. Roman   PetscFunctionBegin;
3392f0a5c5aSJose E. Roman   if (reuse == MAT_INITIAL_MATRIX) {
3402f0a5c5aSJose E. Roman     PetscCall(MatConvert(A, MATAIJ, reuse, B));
3412f0a5c5aSJose E. Roman     PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B));
3422f0a5c5aSJose E. Roman   } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3442f0a5c5aSJose E. Roman }
3452f0a5c5aSJose E. Roman #endif
3462f0a5c5aSJose E. Roman 
3472ef1f0ffSBarry Smith /*MC
3482ef1f0ffSBarry Smith   MATNORMALHERMITIAN - a matrix that behaves like (A*)'*A for `MatMult()` while only containing A
3492ef1f0ffSBarry Smith 
3502ef1f0ffSBarry Smith   Level: intermediate
3512ef1f0ffSBarry Smith 
352*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateNormalHermitian()`, `MatMult()`, `MatNormalHermitianGetMat()`, `MATNORMAL`, `MatCreateNormal()`
3532ef1f0ffSBarry Smith M*/
3542ef1f0ffSBarry Smith 
355ebda224cSfranklin5 /*@
35611a5261eSBarry Smith       MatCreateNormalHermitian - Creates a new matrix object `MATNORMALHERMITIAN` that behaves like (A*)'*A.
357ebda224cSfranklin5 
358c3339decSBarry Smith    Collective
359ebda224cSfranklin5 
360ebda224cSfranklin5    Input Parameter:
361ebda224cSfranklin5 .   A  - the (possibly rectangular complex) matrix
362ebda224cSfranklin5 
363ebda224cSfranklin5    Output Parameter:
364ebda224cSfranklin5 .   N - the matrix that represents (A*)'*A
365ebda224cSfranklin5 
366ebda224cSfranklin5    Level: intermediate
367ebda224cSfranklin5 
36811a5261eSBarry Smith    Note:
36995452b02SPatrick Sanan     The product (A*)'*A is NOT actually formed! Rather the new matrix
37011a5261eSBarry Smith           object performs the matrix-vector product, `MatMult()`, by first multiplying by
371ebda224cSfranklin5           A and then (A*)'
37211a5261eSBarry Smith 
373*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatNormalHermitianGetMat()`
374ebda224cSfranklin5 @*/
375d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormalHermitian(Mat A, Mat *N)
376d71ae5a4SJacob Faibussowitsch {
377ebda224cSfranklin5   PetscInt    m, n;
378ebda224cSfranklin5   Mat_Normal *Na;
3795cb049f5SJose E. Roman   VecType     vtype;
380ebda224cSfranklin5 
381ebda224cSfranklin5   PetscFunctionBegin;
3829566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
3839566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N));
3849566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*N, n, n, PETSC_DECIDE, PETSC_DECIDE));
3859566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMALHERMITIAN));
3869566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap));
3879566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap));
388ebda224cSfranklin5 
3894dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&Na));
390ebda224cSfranklin5   (*N)->data = (void *)Na;
3919566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
392ebda224cSfranklin5   Na->A     = A;
393ebda224cSfranklin5   Na->scale = 1.0;
394ebda224cSfranklin5 
3959566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, NULL, &Na->w));
396ebda224cSfranklin5 
3972e5cc420SJose E. Roman   (*N)->ops->destroy           = MatDestroy_NormalHermitian;
3982e5cc420SJose E. Roman   (*N)->ops->mult              = MatMult_NormalHermitian;
399ebda224cSfranklin5   (*N)->ops->multtranspose     = MatMultHermitianTranspose_Normal;
400ebda224cSfranklin5   (*N)->ops->multtransposeadd  = MatMultHermitianTransposeAdd_Normal;
401ebda224cSfranklin5   (*N)->ops->multadd           = MatMultHermitianAdd_Normal;
4022e5cc420SJose E. Roman   (*N)->ops->getdiagonal       = MatGetDiagonal_NormalHermitian;
4032e5cc420SJose E. Roman   (*N)->ops->getdiagonalblock  = MatGetDiagonalBlock_NormalHermitian;
4042e5cc420SJose E. Roman   (*N)->ops->scale             = MatScale_NormalHermitian;
4052e5cc420SJose E. Roman   (*N)->ops->diagonalscale     = MatDiagonalScale_NormalHermitian;
406ab704866SPierre Jolivet   (*N)->ops->createsubmatrices = MatCreateSubMatrices_NormalHermitian;
407ab704866SPierre Jolivet   (*N)->ops->permute           = MatPermute_NormalHermitian;
408ab704866SPierre Jolivet   (*N)->ops->duplicate         = MatDuplicate_NormalHermitian;
409ab704866SPierre Jolivet   (*N)->ops->copy              = MatCopy_NormalHermitian;
410ebda224cSfranklin5   (*N)->assembled              = PETSC_TRUE;
4115cb049f5SJose E. Roman   (*N)->preallocated           = PETSC_TRUE;
41212ab1b64SPierre Jolivet 
4132e5cc420SJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatNormalGetMatHermitian_C", MatNormalGetMat_NormalHermitian));
414ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_seqaij_C", MatConvert_NormalHermitian_AIJ));
415ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_mpiaij_C", MatConvert_NormalHermitian_AIJ));
4162f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
4172f0a5c5aSJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_hypre_C", MatConvert_NormalHermitian_HYPRE));
4182f0a5c5aSJose E. Roman #endif
4199566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*N, MAT_HERMITIAN, PETSC_TRUE));
4209566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A, &vtype));
4219566063dSJacob Faibussowitsch   PetscCall(MatSetVecType(*N, vtype));
4225cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE)
4239566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(*N, A->boundtocpu));
4245cb049f5SJose E. Roman #endif
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
426ebda224cSfranklin5 }
427