xref: /petsc/src/mat/impls/normal/normmh.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 
119371c9d4SSatish Balay PetscErrorCode MatScale_NormalHermitian(Mat inA, PetscScalar scale) {
12ebda224cSfranklin5   Mat_Normal *a = (Mat_Normal *)inA->data;
13ebda224cSfranklin5 
14ebda224cSfranklin5   PetscFunctionBegin;
15ebda224cSfranklin5   a->scale *= scale;
16ebda224cSfranklin5   PetscFunctionReturn(0);
17ebda224cSfranklin5 }
18ebda224cSfranklin5 
199371c9d4SSatish Balay PetscErrorCode MatDiagonalScale_NormalHermitian(Mat inA, Vec left, Vec right) {
20ebda224cSfranklin5   Mat_Normal *a = (Mat_Normal *)inA->data;
21ebda224cSfranklin5 
22ebda224cSfranklin5   PetscFunctionBegin;
23ebda224cSfranklin5   if (left) {
24ebda224cSfranklin5     if (!a->left) {
259566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(left, &a->left));
269566063dSJacob Faibussowitsch       PetscCall(VecCopy(left, a->left));
27ebda224cSfranklin5     } else {
289566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->left, left, a->left));
29ebda224cSfranklin5     }
30ebda224cSfranklin5   }
31ebda224cSfranklin5   if (right) {
32ebda224cSfranklin5     if (!a->right) {
339566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(right, &a->right));
349566063dSJacob Faibussowitsch       PetscCall(VecCopy(right, a->right));
35ebda224cSfranklin5     } else {
369566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(a->right, right, a->right));
37ebda224cSfranklin5     }
38ebda224cSfranklin5   }
39ebda224cSfranklin5   PetscFunctionReturn(0);
40ebda224cSfranklin5 }
41ebda224cSfranklin5 
429371c9d4SSatish Balay PetscErrorCode MatCreateSubMatrices_NormalHermitian(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) {
43ab704866SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)mat->data;
44ab704866SPierre Jolivet   Mat         B = a->A, *suba;
45ab704866SPierre Jolivet   IS         *row;
46ab704866SPierre Jolivet   PetscInt    M;
47ab704866SPierre Jolivet 
48ab704866SPierre Jolivet   PetscFunctionBegin;
49ab704866SPierre Jolivet   PetscCheck(!a->left && !a->right && irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented");
5048a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
51ab704866SPierre Jolivet   PetscCall(MatGetSize(B, &M, NULL));
52ab704866SPierre Jolivet   PetscCall(PetscMalloc1(n, &row));
53ab704866SPierre Jolivet   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0]));
54ab704866SPierre Jolivet   PetscCall(ISSetIdentity(row[0]));
55ab704866SPierre Jolivet   for (M = 1; M < n; ++M) row[M] = row[0];
56ab704866SPierre Jolivet   PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba));
57ab704866SPierre Jolivet   for (M = 0; M < n; ++M) {
58ab704866SPierre Jolivet     PetscCall(MatCreateNormalHermitian(suba[M], *submat + M));
59ab704866SPierre Jolivet     ((Mat_Normal *)(*submat)[M]->data)->scale = a->scale;
60ab704866SPierre Jolivet   }
61ab704866SPierre Jolivet   PetscCall(ISDestroy(&row[0]));
62ab704866SPierre Jolivet   PetscCall(PetscFree(row));
63ab704866SPierre Jolivet   PetscCall(MatDestroySubMatrices(n, &suba));
64ab704866SPierre Jolivet   PetscFunctionReturn(0);
65ab704866SPierre Jolivet }
66ab704866SPierre Jolivet 
679371c9d4SSatish Balay PetscErrorCode MatPermute_NormalHermitian(Mat A, IS rowp, IS colp, Mat *B) {
68ab704866SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
69ab704866SPierre Jolivet   Mat         C, Aa = a->A;
70ab704866SPierre Jolivet   IS          row;
71ab704866SPierre Jolivet 
72ab704866SPierre Jolivet   PetscFunctionBegin;
73ab704866SPierre Jolivet   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
74ab704866SPierre Jolivet   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row));
75ab704866SPierre Jolivet   PetscCall(ISSetIdentity(row));
76ab704866SPierre Jolivet   PetscCall(MatPermute(Aa, row, colp, &C));
77ab704866SPierre Jolivet   PetscCall(ISDestroy(&row));
78ab704866SPierre Jolivet   PetscCall(MatCreateNormalHermitian(C, B));
79ab704866SPierre Jolivet   PetscCall(MatDestroy(&C));
80ab704866SPierre Jolivet   PetscFunctionReturn(0);
81ab704866SPierre Jolivet }
82ab704866SPierre Jolivet 
839371c9d4SSatish Balay PetscErrorCode MatDuplicate_NormalHermitian(Mat A, MatDuplicateOption op, Mat *B) {
84ab704866SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data;
85ab704866SPierre Jolivet   Mat         C;
86ab704866SPierre Jolivet 
87ab704866SPierre Jolivet   PetscFunctionBegin;
88ab704866SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
89ab704866SPierre Jolivet   PetscCall(MatDuplicate(a->A, op, &C));
90ab704866SPierre Jolivet   PetscCall(MatCreateNormalHermitian(C, B));
91ab704866SPierre Jolivet   PetscCall(MatDestroy(&C));
92ab704866SPierre Jolivet   if (op == MAT_COPY_VALUES) ((Mat_Normal *)(*B)->data)->scale = a->scale;
93ab704866SPierre Jolivet   PetscFunctionReturn(0);
94ab704866SPierre Jolivet }
95ab704866SPierre Jolivet 
969371c9d4SSatish Balay PetscErrorCode MatCopy_NormalHermitian(Mat A, Mat B, MatStructure str) {
97ab704866SPierre Jolivet   Mat_Normal *a = (Mat_Normal *)A->data, *b = (Mat_Normal *)B->data;
98ab704866SPierre Jolivet 
99ab704866SPierre Jolivet   PetscFunctionBegin;
100ab704866SPierre Jolivet   PetscCheck(!a->left && !a->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not implemented");
101ab704866SPierre Jolivet   PetscCall(MatCopy(a->A, b->A, str));
102ab704866SPierre Jolivet   b->scale = a->scale;
103ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->left));
104ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->right));
105ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->leftwork));
106ab704866SPierre Jolivet   PetscCall(VecDestroy(&b->rightwork));
107ab704866SPierre Jolivet   PetscFunctionReturn(0);
108ab704866SPierre Jolivet }
109ab704866SPierre Jolivet 
1109371c9d4SSatish Balay PetscErrorCode MatMult_NormalHermitian(Mat N, Vec x, Vec y) {
111ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
112ebda224cSfranklin5   Vec         in;
113ebda224cSfranklin5 
114ebda224cSfranklin5   PetscFunctionBegin;
115ebda224cSfranklin5   in = x;
116ebda224cSfranklin5   if (Na->right) {
11748a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1189566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
119ebda224cSfranklin5     in = Na->rightwork;
120ebda224cSfranklin5   }
1219566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1229566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y));
1231baa6e33SBarry Smith   if (Na->left) PetscCall(VecPointwiseMult(y, Na->left, y));
1249566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
125ebda224cSfranklin5   PetscFunctionReturn(0);
126ebda224cSfranklin5 }
127ebda224cSfranklin5 
1289371c9d4SSatish Balay PetscErrorCode MatMultHermitianAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3) {
129ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
130ebda224cSfranklin5   Vec         in;
131ebda224cSfranklin5 
132ebda224cSfranklin5   PetscFunctionBegin;
133ebda224cSfranklin5   in = v1;
134ebda224cSfranklin5   if (Na->right) {
13548a46eb9SPierre Jolivet     if (!Na->rightwork) PetscCall(VecDuplicate(Na->right, &Na->rightwork));
1369566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->rightwork, Na->right, in));
137ebda224cSfranklin5     in = Na->rightwork;
138ebda224cSfranklin5   }
1399566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1409566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
141ebda224cSfranklin5   if (Na->left) {
1429566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(Na->A, Na->w, v3));
1439566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->left, v3));
1449566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
145ebda224cSfranklin5   } else {
1469566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTransposeAdd(Na->A, Na->w, v2, v3));
147ebda224cSfranklin5   }
148ebda224cSfranklin5   PetscFunctionReturn(0);
149ebda224cSfranklin5 }
150ebda224cSfranklin5 
1519371c9d4SSatish Balay PetscErrorCode MatMultHermitianTranspose_Normal(Mat N, Vec x, Vec y) {
152ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
153ebda224cSfranklin5   Vec         in;
154ebda224cSfranklin5 
155ebda224cSfranklin5   PetscFunctionBegin;
156ebda224cSfranklin5   in = x;
157ebda224cSfranklin5   if (Na->left) {
15848a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
1599566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
160ebda224cSfranklin5     in = Na->leftwork;
161ebda224cSfranklin5   }
1629566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1639566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y));
1641baa6e33SBarry Smith   if (Na->right) PetscCall(VecPointwiseMult(y, Na->right, y));
1659566063dSJacob Faibussowitsch   PetscCall(VecScale(y, Na->scale));
166ebda224cSfranklin5   PetscFunctionReturn(0);
167ebda224cSfranklin5 }
168ebda224cSfranklin5 
1699371c9d4SSatish Balay PetscErrorCode MatMultHermitianTransposeAdd_Normal(Mat N, Vec v1, Vec v2, Vec v3) {
170ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
171ebda224cSfranklin5   Vec         in;
172ebda224cSfranklin5 
173ebda224cSfranklin5   PetscFunctionBegin;
174ebda224cSfranklin5   in = v1;
175ebda224cSfranklin5   if (Na->left) {
17648a46eb9SPierre Jolivet     if (!Na->leftwork) PetscCall(VecDuplicate(Na->left, &Na->leftwork));
1779566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Na->leftwork, Na->left, in));
178ebda224cSfranklin5     in = Na->leftwork;
179ebda224cSfranklin5   }
1809566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A, in, Na->w));
1819566063dSJacob Faibussowitsch   PetscCall(VecScale(Na->w, Na->scale));
182ebda224cSfranklin5   if (Na->right) {
1839566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTranspose(Na->A, Na->w, v3));
1849566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(v3, Na->right, v3));
1859566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v3, 1.0, v2));
186ebda224cSfranklin5   } else {
1879566063dSJacob Faibussowitsch     PetscCall(MatMultHermitianTransposeAdd(Na->A, Na->w, v2, v3));
188ebda224cSfranklin5   }
189ebda224cSfranklin5   PetscFunctionReturn(0);
190ebda224cSfranklin5 }
191ebda224cSfranklin5 
1929371c9d4SSatish Balay PetscErrorCode MatDestroy_NormalHermitian(Mat N) {
193ebda224cSfranklin5   Mat_Normal *Na = (Mat_Normal *)N->data;
194ebda224cSfranklin5 
195ebda224cSfranklin5   PetscFunctionBegin;
1969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
197a48507d3SJose E. Roman   PetscCall(MatDestroy(&Na->D));
1989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->w));
1999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->left));
2009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->right));
2019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->leftwork));
2029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->rightwork));
2039566063dSJacob Faibussowitsch   PetscCall(PetscFree(N->data));
2049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMatHermitian_C", NULL));
205ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_seqaij_C", NULL));
206ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_mpiaij_C", NULL));
207ebda224cSfranklin5   PetscFunctionReturn(0);
208ebda224cSfranklin5 }
209ebda224cSfranklin5 
210ebda224cSfranklin5 /*
211ebda224cSfranklin5       Slow, nonscalable version
212ebda224cSfranklin5 */
2139371c9d4SSatish Balay PetscErrorCode MatGetDiagonal_NormalHermitian(Mat N, Vec v) {
214ebda224cSfranklin5   Mat_Normal        *Na = (Mat_Normal *)N->data;
215ebda224cSfranklin5   Mat                A  = Na->A;
216ebda224cSfranklin5   PetscInt           i, j, rstart, rend, nnz;
217ebda224cSfranklin5   const PetscInt    *cols;
218ebda224cSfranklin5   PetscScalar       *diag, *work, *values;
219ebda224cSfranklin5   const PetscScalar *mvalues;
220ebda224cSfranklin5 
221ebda224cSfranklin5   PetscFunctionBegin;
2229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work));
2239566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(work, A->cmap->N));
2249566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
225ebda224cSfranklin5   for (i = rstart; i < rend; i++) {
2269566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues));
227ad540459SPierre Jolivet     for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * PetscConj(mvalues[j]);
2289566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues));
229ebda224cSfranklin5   }
2301c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N)));
231ebda224cSfranklin5   rstart = N->cmap->rstart;
232ebda224cSfranklin5   rend   = N->cmap->rend;
2339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &values));
2349566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart));
2359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &values));
2369566063dSJacob Faibussowitsch   PetscCall(PetscFree2(diag, work));
2379566063dSJacob Faibussowitsch   PetscCall(VecScale(v, Na->scale));
238ebda224cSfranklin5   PetscFunctionReturn(0);
239ebda224cSfranklin5 }
240ebda224cSfranklin5 
2419371c9d4SSatish Balay PetscErrorCode MatGetDiagonalBlock_NormalHermitian(Mat N, Mat *D) {
242a48507d3SJose E. Roman   Mat_Normal *Na = (Mat_Normal *)N->data;
243a48507d3SJose E. Roman   Mat         M, A = Na->A;
244a48507d3SJose E. Roman 
245a48507d3SJose E. Roman   PetscFunctionBegin;
246a48507d3SJose E. Roman   PetscCall(MatGetDiagonalBlock(A, &M));
247a48507d3SJose E. Roman   PetscCall(MatCreateNormalHermitian(M, &Na->D));
248a48507d3SJose E. Roman   *D = Na->D;
249a48507d3SJose E. Roman   PetscFunctionReturn(0);
250a48507d3SJose E. Roman }
251a48507d3SJose E. Roman 
2529371c9d4SSatish Balay PetscErrorCode MatNormalGetMat_NormalHermitian(Mat A, Mat *M) {
25365f45395SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal *)A->data;
25465f45395SPierre Jolivet 
25565f45395SPierre Jolivet   PetscFunctionBegin;
25665f45395SPierre Jolivet   *M = Aa->A;
25765f45395SPierre Jolivet   PetscFunctionReturn(0);
25865f45395SPierre Jolivet }
25965f45395SPierre Jolivet 
26065f45395SPierre Jolivet /*@
26111a5261eSBarry Smith       MatNormalHermitianGetMat - Gets the `Mat` object stored inside a `MATNORMALHERMITIAN`
26265f45395SPierre Jolivet 
26311a5261eSBarry Smith    Logically collective on A
26465f45395SPierre Jolivet 
26565f45395SPierre Jolivet    Input Parameter:
26611a5261eSBarry Smith .   A  - the `MATNORMALHERMITIAN` matrix
26765f45395SPierre Jolivet 
26865f45395SPierre Jolivet    Output Parameter:
26965f45395SPierre Jolivet .   M - the matrix object stored inside A
27065f45395SPierre Jolivet 
27165f45395SPierre Jolivet    Level: intermediate
27265f45395SPierre Jolivet 
27311a5261eSBarry Smith .seealso: `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
27465f45395SPierre Jolivet @*/
2759371c9d4SSatish Balay PetscErrorCode MatNormalHermitianGetMat(Mat A, Mat *M) {
27665f45395SPierre Jolivet   PetscFunctionBegin;
27765f45395SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
27865f45395SPierre Jolivet   PetscValidType(A, 1);
27965f45395SPierre Jolivet   PetscValidPointer(M, 2);
280cac4c232SBarry Smith   PetscUseMethod(A, "MatNormalGetMatHermitian_C", (Mat, Mat *), (A, M));
28165f45395SPierre Jolivet   PetscFunctionReturn(0);
28265f45395SPierre Jolivet }
28365f45395SPierre Jolivet 
2849371c9d4SSatish Balay PetscErrorCode MatConvert_NormalHermitian_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) {
285ab704866SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal *)A->data;
286ab704866SPierre Jolivet   Mat         B, conjugate;
287ab704866SPierre Jolivet   PetscInt    m, n, M, N;
288ab704866SPierre Jolivet 
289ab704866SPierre Jolivet   PetscFunctionBegin;
290ab704866SPierre Jolivet   PetscCall(MatGetSize(A, &M, &N));
291ab704866SPierre Jolivet   PetscCall(MatGetLocalSize(A, &m, &n));
292ab704866SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
293ab704866SPierre Jolivet     B = *newmat;
294ab704866SPierre Jolivet     PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B));
295ab704866SPierre Jolivet   } else {
296ab704866SPierre Jolivet     PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B));
297ab704866SPierre Jolivet     PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
298ab704866SPierre Jolivet     PetscCall(MatProductSetFromOptions(B));
299ab704866SPierre Jolivet     PetscCall(MatProductSymbolic(B));
300ab704866SPierre Jolivet     PetscCall(MatSetOption(B, !PetscDefined(USE_COMPLEX) ? MAT_SYMMETRIC : MAT_HERMITIAN, PETSC_TRUE));
301ab704866SPierre Jolivet   }
302ab704866SPierre Jolivet   if (PetscDefined(USE_COMPLEX)) {
303ab704866SPierre Jolivet     PetscCall(MatDuplicate(Aa->A, MAT_COPY_VALUES, &conjugate));
304ab704866SPierre Jolivet     PetscCall(MatConjugate(conjugate));
305ab704866SPierre Jolivet     PetscCall(MatProductReplaceMats(conjugate, Aa->A, NULL, B));
306ab704866SPierre Jolivet   }
307ab704866SPierre Jolivet   PetscCall(MatProductNumeric(B));
308ab704866SPierre Jolivet   if (PetscDefined(USE_COMPLEX)) PetscCall(MatDestroy(&conjugate));
309ab704866SPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
310ab704866SPierre Jolivet     PetscCall(MatHeaderReplace(A, &B));
311ab704866SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
312ab704866SPierre Jolivet   PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat));
313ab704866SPierre Jolivet   PetscFunctionReturn(0);
314ab704866SPierre Jolivet }
315ab704866SPierre Jolivet 
316ebda224cSfranklin5 /*@
31711a5261eSBarry Smith       MatCreateNormalHermitian - Creates a new matrix object `MATNORMALHERMITIAN` that behaves like (A*)'*A.
318ebda224cSfranklin5 
31911a5261eSBarry Smith    Collective on A
320ebda224cSfranklin5 
321ebda224cSfranklin5    Input Parameter:
322ebda224cSfranklin5 .   A  - the (possibly rectangular complex) matrix
323ebda224cSfranklin5 
324ebda224cSfranklin5    Output Parameter:
325ebda224cSfranklin5 .   N - the matrix that represents (A*)'*A
326ebda224cSfranklin5 
327ebda224cSfranklin5    Level: intermediate
328ebda224cSfranklin5 
32911a5261eSBarry Smith    Note:
33095452b02SPatrick Sanan     The product (A*)'*A is NOT actually formed! Rather the new matrix
33111a5261eSBarry Smith           object performs the matrix-vector product, `MatMult()`, by first multiplying by
332ebda224cSfranklin5           A and then (A*)'
33311a5261eSBarry Smith 
33411a5261eSBarry Smith .seealso: `MATNORMAL`, `MATNORMALHERMITIAN`, `MatNormalHermitianGetMat()`
335ebda224cSfranklin5 @*/
3369371c9d4SSatish Balay PetscErrorCode MatCreateNormalHermitian(Mat A, Mat *N) {
337ebda224cSfranklin5   PetscInt    m, n;
338ebda224cSfranklin5   Mat_Normal *Na;
3395cb049f5SJose E. Roman   VecType     vtype;
340ebda224cSfranklin5 
341ebda224cSfranklin5   PetscFunctionBegin;
3429566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
3439566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N));
3449566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*N, n, n, PETSC_DECIDE, PETSC_DECIDE));
3459566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMALHERMITIAN));
3469566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap));
3479566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap));
348ebda224cSfranklin5 
349*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&Na));
350ebda224cSfranklin5   (*N)->data = (void *)Na;
3519566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
352ebda224cSfranklin5   Na->A     = A;
353ebda224cSfranklin5   Na->scale = 1.0;
354ebda224cSfranklin5 
3559566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, NULL, &Na->w));
356ebda224cSfranklin5 
3572e5cc420SJose E. Roman   (*N)->ops->destroy           = MatDestroy_NormalHermitian;
3582e5cc420SJose E. Roman   (*N)->ops->mult              = MatMult_NormalHermitian;
359ebda224cSfranklin5   (*N)->ops->multtranspose     = MatMultHermitianTranspose_Normal;
360ebda224cSfranklin5   (*N)->ops->multtransposeadd  = MatMultHermitianTransposeAdd_Normal;
361ebda224cSfranklin5   (*N)->ops->multadd           = MatMultHermitianAdd_Normal;
3622e5cc420SJose E. Roman   (*N)->ops->getdiagonal       = MatGetDiagonal_NormalHermitian;
3632e5cc420SJose E. Roman   (*N)->ops->getdiagonalblock  = MatGetDiagonalBlock_NormalHermitian;
3642e5cc420SJose E. Roman   (*N)->ops->scale             = MatScale_NormalHermitian;
3652e5cc420SJose E. Roman   (*N)->ops->diagonalscale     = MatDiagonalScale_NormalHermitian;
366ab704866SPierre Jolivet   (*N)->ops->createsubmatrices = MatCreateSubMatrices_NormalHermitian;
367ab704866SPierre Jolivet   (*N)->ops->permute           = MatPermute_NormalHermitian;
368ab704866SPierre Jolivet   (*N)->ops->duplicate         = MatDuplicate_NormalHermitian;
369ab704866SPierre Jolivet   (*N)->ops->copy              = MatCopy_NormalHermitian;
370ebda224cSfranklin5   (*N)->assembled              = PETSC_TRUE;
3715cb049f5SJose E. Roman   (*N)->preallocated           = PETSC_TRUE;
37212ab1b64SPierre Jolivet 
3732e5cc420SJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatNormalGetMatHermitian_C", MatNormalGetMat_NormalHermitian));
374ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_seqaij_C", MatConvert_NormalHermitian_AIJ));
375ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)(*N), "MatConvert_normalh_mpiaij_C", MatConvert_NormalHermitian_AIJ));
3769566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*N, MAT_HERMITIAN, PETSC_TRUE));
3779566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A, &vtype));
3789566063dSJacob Faibussowitsch   PetscCall(MatSetVecType(*N, vtype));
3795cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE)
3809566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(*N, A->boundtocpu));
3815cb049f5SJose E. Roman #endif
382ebda224cSfranklin5   PetscFunctionReturn(0);
383ebda224cSfranklin5 }
384