xref: /petsc/src/mat/impls/normal/normmh.c (revision 14e4dea26ebf7cc95ab4a0045be33d3dec6b4ea6)
1ebda224cSfranklin5 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
2ebda224cSfranklin5 
3ebda224cSfranklin5 typedef struct {
4ebda224cSfranklin5   Mat         A;
5a48507d3SJose E. Roman   Mat         D; /* local submatrix for diagonal part */
6ebda224cSfranklin5   Vec         w, left, right, leftwork, rightwork;
7ebda224cSfranklin5   PetscScalar scale;
813b343c1SBarry Smith } Mat_NormalHermitian;
9ebda224cSfranklin5 
1066976f2fSJacob Faibussowitsch static PetscErrorCode MatScale_NormalHermitian(Mat inA, PetscScalar scale)
11d71ae5a4SJacob Faibussowitsch {
1213b343c1SBarry Smith   Mat_NormalHermitian *a = (Mat_NormalHermitian *)inA->data;
13ebda224cSfranklin5 
14ebda224cSfranklin5   PetscFunctionBegin;
15ebda224cSfranklin5   a->scale *= scale;
163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17ebda224cSfranklin5 }
18ebda224cSfranklin5 
1966976f2fSJacob Faibussowitsch static PetscErrorCode MatDiagonalScale_NormalHermitian(Mat inA, Vec left, Vec right)
20d71ae5a4SJacob Faibussowitsch {
2113b343c1SBarry Smith   Mat_NormalHermitian *a = (Mat_NormalHermitian *)inA->data;
22ebda224cSfranklin5 
23ebda224cSfranklin5   PetscFunctionBegin;
24ebda224cSfranklin5   if (left) {
25ebda224cSfranklin5     if (!a->left) {
269566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(left, &a->left));
279566063dSJacob Faibussowitsch       PetscCall(VecCopy(left, a->left));
28ebda224cSfranklin5     } else {
299566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->left, left, a->left));
30ebda224cSfranklin5     }
31ebda224cSfranklin5   }
32ebda224cSfranklin5   if (right) {
33ebda224cSfranklin5     if (!a->right) {
349566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(right, &a->right));
359566063dSJacob Faibussowitsch       PetscCall(VecCopy(right, a->right));
36ebda224cSfranklin5     } else {
379566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->right, right, a->right));
38ebda224cSfranklin5     }
39ebda224cSfranklin5   }
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41ebda224cSfranklin5 }
42ebda224cSfranklin5 
4366976f2fSJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_NormalHermitian(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
44d71ae5a4SJacob Faibussowitsch {
4513b343c1SBarry Smith   Mat_NormalHermitian *a = (Mat_NormalHermitian *)mat->data;
46ab704866SPierre Jolivet   Mat                  B = a->A, *suba;
47ab704866SPierre Jolivet   IS                  *row;
48ab704866SPierre Jolivet   PetscInt             M;
49ab704866SPierre Jolivet 
50ab704866SPierre Jolivet   PetscFunctionBegin;
51ab704866SPierre Jolivet   PetscCheck(!a->left && !a->right && irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented");
5248a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
53ab704866SPierre Jolivet   PetscCall(MatGetSize(B, &M, NULL));
54ab704866SPierre Jolivet   PetscCall(PetscMalloc1(n, &row));
55ab704866SPierre Jolivet   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0]));
56ab704866SPierre Jolivet   PetscCall(ISSetIdentity(row[0]));
57ab704866SPierre Jolivet   for (M = 1; M < n; ++M) row[M] = row[0];
58ab704866SPierre Jolivet   PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba));
59ab704866SPierre Jolivet   for (M = 0; M < n; ++M) {
60ab704866SPierre Jolivet     PetscCall(MatCreateNormalHermitian(suba[M], *submat + M));
6113b343c1SBarry Smith     ((Mat_NormalHermitian *)(*submat)[M]->data)->scale = a->scale;
62ab704866SPierre Jolivet   }
63ab704866SPierre Jolivet   PetscCall(ISDestroy(&row[0]));
64ab704866SPierre Jolivet   PetscCall(PetscFree(row));
65ab704866SPierre Jolivet   PetscCall(MatDestroySubMatrices(n, &suba));
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67ab704866SPierre Jolivet }
68ab704866SPierre Jolivet 
6966976f2fSJacob Faibussowitsch static PetscErrorCode MatPermute_NormalHermitian(Mat A, IS rowp, IS colp, Mat *B)
70d71ae5a4SJacob Faibussowitsch {
7113b343c1SBarry Smith   Mat_NormalHermitian *a = (Mat_NormalHermitian *)A->data;
72ab704866SPierre Jolivet   Mat                  C, Aa = a->A;
73ab704866SPierre Jolivet   IS                   row;
74ab704866SPierre Jolivet 
75ab704866SPierre Jolivet   PetscFunctionBegin;
76ab704866SPierre Jolivet   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
77ab704866SPierre Jolivet   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row));
78ab704866SPierre Jolivet   PetscCall(ISSetIdentity(row));
79ab704866SPierre Jolivet   PetscCall(MatPermute(Aa, row, colp, &C));
80ab704866SPierre Jolivet   PetscCall(ISDestroy(&row));
81ab704866SPierre Jolivet   PetscCall(MatCreateNormalHermitian(C, B));
82ab704866SPierre Jolivet   PetscCall(MatDestroy(&C));
833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84ab704866SPierre Jolivet }
85ab704866SPierre Jolivet 
8666976f2fSJacob Faibussowitsch static PetscErrorCode MatDuplicate_NormalHermitian(Mat A, MatDuplicateOption op, Mat *B)
87d71ae5a4SJacob Faibussowitsch {
8813b343c1SBarry Smith   Mat_NormalHermitian *a = (Mat_NormalHermitian *)A->data;
89ab704866SPierre Jolivet   Mat                  C;
90ab704866SPierre Jolivet 
91ab704866SPierre Jolivet   PetscFunctionBegin;
92ab704866SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
93ab704866SPierre Jolivet   PetscCall(MatDuplicate(a->A, op, &C));
94ab704866SPierre Jolivet   PetscCall(MatCreateNormalHermitian(C, B));
95ab704866SPierre Jolivet   PetscCall(MatDestroy(&C));
9613b343c1SBarry Smith   if (op == MAT_COPY_VALUES) ((Mat_NormalHermitian *)(*B)->data)->scale = a->scale;
973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98ab704866SPierre Jolivet }
99ab704866SPierre Jolivet 
10066976f2fSJacob Faibussowitsch static PetscErrorCode MatCopy_NormalHermitian(Mat A, Mat B, MatStructure str)
101d71ae5a4SJacob Faibussowitsch {
10213b343c1SBarry Smith   Mat_NormalHermitian *a = (Mat_NormalHermitian *)A->data, *b = (Mat_NormalHermitian *)B->data;
103ab704866SPierre Jolivet 
104ab704866SPierre Jolivet   PetscFunctionBegin;
105ab704866SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
106ab704866SPierre Jolivet   PetscCall(MatCopy(a->A, b->A, str));
107ab704866SPierre Jolivet   b->scale = a->scale;
108ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->left));
109ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->right));
110ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->leftwork));
111ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->rightwork));
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113ab704866SPierre Jolivet }
114ab704866SPierre Jolivet 
11566976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_NormalHermitian(Mat N, Vec x, Vec y)
116d71ae5a4SJacob Faibussowitsch {
11713b343c1SBarry Smith   Mat_NormalHermitian *Na = (Mat_NormalHermitian *)N->data;
118ebda224cSfranklin5   Vec                  in;
119ebda224cSfranklin5 
120ebda224cSfranklin5   PetscFunctionBegin;
121ebda224cSfranklin5   in = x;
122ebda224cSfranklin5   if (Na->right) {
12348a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1249566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
125ebda224cSfranklin5     in = Na->rightwork;
126ebda224cSfranklin5   }
1279566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1289566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y));
1291baa6e33SBarry Smith   if (Na->left) PetscCall(VecPointwiseMult(y, Na->left, y));
1309566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
1313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
132ebda224cSfranklin5 }
133ebda224cSfranklin5 
13466976f2fSJacob Faibussowitsch static PetscErrorCode MatMultHermitianAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3)
135d71ae5a4SJacob Faibussowitsch {
13613b343c1SBarry Smith   Mat_NormalHermitian *Na = (Mat_NormalHermitian *)N->data;
137ebda224cSfranklin5   Vec                  in;
138ebda224cSfranklin5 
139ebda224cSfranklin5   PetscFunctionBegin;
140ebda224cSfranklin5   in = v1;
141ebda224cSfranklin5   if (Na->right) {
14248a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1439566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
144ebda224cSfranklin5     in = Na->rightwork;
145ebda224cSfranklin5   }
1469566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1479566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
148ebda224cSfranklin5   if (Na->left) {
1499566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(Na->A, Na->w, v3));
1509566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->left, v3));
1519566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
152ebda224cSfranklin5   } else {
1539566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTransposeAdd(Na->A, Na->w, v2, v3));
154ebda224cSfranklin5   }
1553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
156ebda224cSfranklin5 }
157ebda224cSfranklin5 
15866976f2fSJacob Faibussowitsch static PetscErrorCode MatMultHermitianTranspose_Normal(Mat N, Vec x, Vec y)
159d71ae5a4SJacob Faibussowitsch {
16013b343c1SBarry Smith   Mat_NormalHermitian *Na = (Mat_NormalHermitian *)N->data;
161ebda224cSfranklin5   Vec                  in;
162ebda224cSfranklin5 
163ebda224cSfranklin5   PetscFunctionBegin;
164ebda224cSfranklin5   in = x;
165ebda224cSfranklin5   if (Na->left) {
16648a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
1679566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
168ebda224cSfranklin5     in = Na->leftwork;
169ebda224cSfranklin5   }
1709566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1719566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y));
1721baa6e33SBarry Smith   if (Na->right) PetscCall(VecPointwiseMult(y, Na->right, y));
1739566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
1743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
175ebda224cSfranklin5 }
176ebda224cSfranklin5 
17766976f2fSJacob Faibussowitsch static PetscErrorCode MatMultHermitianTransposeAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3)
178d71ae5a4SJacob Faibussowitsch {
17913b343c1SBarry Smith   Mat_NormalHermitian *Na = (Mat_NormalHermitian *)N->data;
180ebda224cSfranklin5   Vec                  in;
181ebda224cSfranklin5 
182ebda224cSfranklin5   PetscFunctionBegin;
183ebda224cSfranklin5   in = v1;
184ebda224cSfranklin5   if (Na->left) {
18548a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
1869566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
187ebda224cSfranklin5     in = Na->leftwork;
188ebda224cSfranklin5   }
1899566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1909566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
191ebda224cSfranklin5   if (Na->right) {
1929566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(Na->A, Na->w, v3));
1939566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->right, v3));
1949566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
195ebda224cSfranklin5   } else {
1969566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTransposeAdd(Na->A, Na->w, v2, v3));
197ebda224cSfranklin5   }
1983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
199ebda224cSfranklin5 }
200ebda224cSfranklin5 
20166976f2fSJacob Faibussowitsch static PetscErrorCode MatDestroy_NormalHermitian(Mat N)
202d71ae5a4SJacob Faibussowitsch {
20313b343c1SBarry Smith   Mat_NormalHermitian *Na = (Mat_NormalHermitian *)N->data;
204ebda224cSfranklin5 
205ebda224cSfranklin5   PetscFunctionBegin;
2069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
207a48507d3SJose E. Roman   PetscCall(MatDestroy(&Na->D));
2089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->w));
2099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->left));
2109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->right));
2119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->leftwork));
2129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->rightwork));
2139566063dSJacob Faibussowitsch   PetscCall(PetscFree(N->data));
214*14e4dea2SJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalHermitianGetMat_C", NULL));
215ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_seqaij_C", NULL));
216ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_mpiaij_C", NULL));
2172f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
2182f0a5c5aSJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_hypre_C", NULL));
2192f0a5c5aSJose E. Roman #endif
2203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
221ebda224cSfranklin5 }
222ebda224cSfranklin5 
223ebda224cSfranklin5 /*
224ebda224cSfranklin5       Slow, nonscalable version
225ebda224cSfranklin5 */
22666976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_NormalHermitian(Mat N, Vec v)
227d71ae5a4SJacob Faibussowitsch {
22813b343c1SBarry Smith   Mat_NormalHermitian *Na = (Mat_NormalHermitian *)N->data;
229ebda224cSfranklin5   Mat                  A  = Na->A;
230ebda224cSfranklin5   PetscInt             i, j, rstart, rend, nnz;
231ebda224cSfranklin5   const PetscInt      *cols;
232ebda224cSfranklin5   PetscScalar         *diag, *work, *values;
233ebda224cSfranklin5   const PetscScalar   *mvalues;
234ebda224cSfranklin5 
235ebda224cSfranklin5   PetscFunctionBegin;
2369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work));
2379566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(work, A->cmap->N));
2389566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
239ebda224cSfranklin5   for (i = rstart; i < rend; i++) {
2409566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues));
241ad540459SPierre Jolivet     for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * PetscConj(mvalues[j]);
2429566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues));
243ebda224cSfranklin5   }
2441c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N)));
245ebda224cSfranklin5   rstart = N->cmap->rstart;
246ebda224cSfranklin5   rend   = N->cmap->rend;
2479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &values));
2489566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart));
2499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &values));
2509566063dSJacob Faibussowitsch   PetscCall(PetscFree2(diag, work));
2519566063dSJacob Faibussowitsch   PetscCall(VecScale(v, Na->scale));
2523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
253ebda224cSfranklin5 }
254ebda224cSfranklin5 
25566976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_NormalHermitian(Mat N, Mat *D)
256d71ae5a4SJacob Faibussowitsch {
25713b343c1SBarry Smith   Mat_NormalHermitian *Na = (Mat_NormalHermitian *)N->data;
258a48507d3SJose E. Roman   Mat                  M, A = Na->A;
259a48507d3SJose E. Roman 
260a48507d3SJose E. Roman   PetscFunctionBegin;
261a48507d3SJose E. Roman   PetscCall(MatGetDiagonalBlock(A, &M));
262a48507d3SJose E. Roman   PetscCall(MatCreateNormalHermitian(M, &Na->D));
263a48507d3SJose E. Roman   *D = Na->D;
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
265a48507d3SJose E. Roman }
266a48507d3SJose E. Roman 
267*14e4dea2SJose E. Roman static PetscErrorCode MatNormalHermitianGetMat_NormalHermitian(Mat A, Mat *M)
268d71ae5a4SJacob Faibussowitsch {
26913b343c1SBarry Smith   Mat_NormalHermitian *Aa = (Mat_NormalHermitian *)A->data;
27065f45395SPierre Jolivet 
27165f45395SPierre Jolivet   PetscFunctionBegin;
27265f45395SPierre Jolivet   *M = Aa->A;
2733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27465f45395SPierre Jolivet }
27565f45395SPierre Jolivet 
27665f45395SPierre Jolivet /*@
27711a5261eSBarry Smith   MatNormalHermitianGetMat - Gets the `Mat` object stored inside a `MATNORMALHERMITIAN`
27865f45395SPierre Jolivet 
2792ef1f0ffSBarry Smith   Logically Collective
28065f45395SPierre Jolivet 
28165f45395SPierre Jolivet   Input Parameter:
28211a5261eSBarry Smith . A - the `MATNORMALHERMITIAN` matrix
28365f45395SPierre Jolivet 
28465f45395SPierre Jolivet   Output Parameter:
28565f45395SPierre Jolivet . M - the matrix object stored inside A
28665f45395SPierre Jolivet 
28765f45395SPierre Jolivet   Level: intermediate
28865f45395SPierre Jolivet 
2891cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
29065f45395SPierre Jolivet @*/
291d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalHermitianGetMat(Mat A, Mat *M)
292d71ae5a4SJacob Faibussowitsch {
29365f45395SPierre Jolivet   PetscFunctionBegin;
29465f45395SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
29565f45395SPierre Jolivet   PetscValidType(A, 1);
2964f572ea9SToby Isaac   PetscAssertPointer(M, 2);
297*14e4dea2SJose E. Roman   PetscUseMethod(A, "MatNormalHermitianGetMat_C", (Mat, Mat *), (A, M));
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29965f45395SPierre Jolivet }
30065f45395SPierre Jolivet 
30166976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_NormalHermitian_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
302d71ae5a4SJacob Faibussowitsch {
30313b343c1SBarry Smith   Mat_NormalHermitian *Aa = (Mat_NormalHermitian *)A->data;
304ab704866SPierre Jolivet   Mat                  B, conjugate;
305ab704866SPierre Jolivet   PetscInt             m, n, M, N;
306ab704866SPierre Jolivet 
307ab704866SPierre Jolivet   PetscFunctionBegin;
308ab704866SPierre Jolivet   PetscCall(MatGetSize(A, &M, &N));
309ab704866SPierre Jolivet   PetscCall(MatGetLocalSize(A, &m, &n));
310ab704866SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
311ab704866SPierre Jolivet     B = *newmat;
312ab704866SPierre Jolivet     PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B));
313ab704866SPierre Jolivet   } else {
314ab704866SPierre Jolivet     PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B));
315ab704866SPierre Jolivet     PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
316ab704866SPierre Jolivet     PetscCall(MatProductSetFromOptions(B));
317ab704866SPierre Jolivet     PetscCall(MatProductSymbolic(B));
318ab704866SPierre Jolivet     PetscCall(MatSetOption(B, !PetscDefined(USE_COMPLEX) ? MAT_SYMMETRIC : MAT_HERMITIAN, PETSC_TRUE));
319ab704866SPierre Jolivet   }
320ab704866SPierre Jolivet   if (PetscDefined(USE_COMPLEX)) {
321ab704866SPierre Jolivet     PetscCall(MatDuplicate(Aa->A, MAT_COPY_VALUES, &conjugate));
322ab704866SPierre Jolivet     PetscCall(MatConjugate(conjugate));
323ab704866SPierre Jolivet     PetscCall(MatProductReplaceMats(conjugate, Aa->A, NULL, B));
324ab704866SPierre Jolivet   }
325ab704866SPierre Jolivet   PetscCall(MatProductNumeric(B));
326ab704866SPierre Jolivet   if (PetscDefined(USE_COMPLEX)) PetscCall(MatDestroy(&conjugate));
327ab704866SPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
328ab704866SPierre Jolivet     PetscCall(MatHeaderReplace(A, &B));
329ab704866SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
330ab704866SPierre Jolivet   PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat));
3313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
332ab704866SPierre Jolivet }
333ab704866SPierre Jolivet 
3342f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
33566976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_NormalHermitian_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
3362f0a5c5aSJose E. Roman {
3372f0a5c5aSJose E. Roman   PetscFunctionBegin;
3382f0a5c5aSJose E. Roman   if (reuse == MAT_INITIAL_MATRIX) {
3392f0a5c5aSJose E. Roman     PetscCall(MatConvert(A, MATAIJ, reuse, B));
3402f0a5c5aSJose E. Roman     PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B));
3412f0a5c5aSJose E. Roman   } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */
3423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3432f0a5c5aSJose E. Roman }
3442f0a5c5aSJose E. Roman #endif
3452f0a5c5aSJose E. Roman 
3462ef1f0ffSBarry Smith /*MC
3472ef1f0ffSBarry Smith   MATNORMALHERMITIAN - a matrix that behaves like (A*)'*A for `MatMult()` while only containing A
3482ef1f0ffSBarry Smith 
3492ef1f0ffSBarry Smith   Level: intermediate
3502ef1f0ffSBarry Smith 
3511cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateNormalHermitian()`, `MatMult()`, `MatNormalHermitianGetMat()`, `MATNORMAL`, `MatCreateNormal()`
3522ef1f0ffSBarry Smith M*/
3532ef1f0ffSBarry Smith 
354ebda224cSfranklin5 /*@
35511a5261eSBarry Smith   MatCreateNormalHermitian - Creates a new matrix object `MATNORMALHERMITIAN` that behaves like (A*)'*A.
356ebda224cSfranklin5 
357c3339decSBarry Smith   Collective
358ebda224cSfranklin5 
359ebda224cSfranklin5   Input Parameter:
360ebda224cSfranklin5 . A - the (possibly rectangular complex) matrix
361ebda224cSfranklin5 
362ebda224cSfranklin5   Output Parameter:
363ebda224cSfranklin5 . N - the matrix that represents (A*)'*A
364ebda224cSfranklin5 
365ebda224cSfranklin5   Level: intermediate
366ebda224cSfranklin5 
36711a5261eSBarry Smith   Note:
36895452b02SPatrick Sanan   The product (A*)'*A is NOT actually formed! Rather the new matrix
36911a5261eSBarry Smith   object performs the matrix-vector product, `MatMult()`, by first multiplying by
370ebda224cSfranklin5   A and then (A*)'
37111a5261eSBarry Smith 
3721cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatNormalHermitianGetMat()`
373ebda224cSfranklin5 @*/
374d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormalHermitian(Mat A, Mat *N)
375d71ae5a4SJacob Faibussowitsch {
376ebda224cSfranklin5   PetscInt             m, n;
37713b343c1SBarry Smith   Mat_NormalHermitian *Na;
3785cb049f5SJose E. Roman   VecType              vtype;
379ebda224cSfranklin5 
380ebda224cSfranklin5   PetscFunctionBegin;
3819566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
3829566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N));
3839566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*N, n, n, PETSC_DECIDE, PETSC_DECIDE));
3849566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMALHERMITIAN));
3859566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap));
3869566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap));
387ebda224cSfranklin5 
3884dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&Na));
389ebda224cSfranklin5   (*N)->data = (void *)Na;
3909566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
391ebda224cSfranklin5   Na->A     = A;
392ebda224cSfranklin5   Na->scale = 1.0;
393ebda224cSfranklin5 
3949566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, NULL, &Na->w));
395ebda224cSfranklin5 
3962e5cc420SJose E. Roman   (*N)->ops->destroy           = MatDestroy_NormalHermitian;
3972e5cc420SJose E. Roman   (*N)->ops->mult              = MatMult_NormalHermitian;
398ebda224cSfranklin5   (*N)->ops->multtranspose     = MatMultHermitianTranspose_Normal;
399ebda224cSfranklin5   (*N)->ops->multtransposeadd  = MatMultHermitianTransposeAdd_Normal;
400ebda224cSfranklin5   (*N)->ops->multadd           = MatMultHermitianAdd_Normal;
4012e5cc420SJose E. Roman   (*N)->ops->getdiagonal       = MatGetDiagonal_NormalHermitian;
4022e5cc420SJose E. Roman   (*N)->ops->getdiagonalblock  = MatGetDiagonalBlock_NormalHermitian;
4032e5cc420SJose E. Roman   (*N)->ops->scale             = MatScale_NormalHermitian;
4042e5cc420SJose E. Roman   (*N)->ops->diagonalscale     = MatDiagonalScale_NormalHermitian;
405ab704866SPierre Jolivet   (*N)->ops->createsubmatrices = MatCreateSubMatrices_NormalHermitian;
406ab704866SPierre Jolivet   (*N)->ops->permute           = MatPermute_NormalHermitian;
407ab704866SPierre Jolivet   (*N)->ops->duplicate         = MatDuplicate_NormalHermitian;
408ab704866SPierre Jolivet   (*N)->ops->copy              = MatCopy_NormalHermitian;
409ebda224cSfranklin5   (*N)->assembled              = PETSC_TRUE;
4105cb049f5SJose E. Roman   (*N)->preallocated           = PETSC_TRUE;
41112ab1b64SPierre Jolivet 
412*14e4dea2SJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatNormalHermitianGetMat_C", MatNormalHermitianGetMat_NormalHermitian));
413ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_seqaij_C", MatConvert_NormalHermitian_AIJ));
414ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_mpiaij_C", MatConvert_NormalHermitian_AIJ));
4152f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
4162f0a5c5aSJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_hypre_C", MatConvert_NormalHermitian_HYPRE));
4172f0a5c5aSJose E. Roman #endif
4189566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*N, MAT_HERMITIAN, PETSC_TRUE));
4199566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A, &vtype));
4209566063dSJacob Faibussowitsch   PetscCall(MatSetVecType(*N, vtype));
4215cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE)
4229566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(*N, A->boundtocpu));
4235cb049f5SJose E. Roman #endif
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
425ebda224cSfranklin5 }
426