xref: /petsc/src/mat/impls/normal/normmh.c (revision b9c875b8ccb0466cfbdbcc08c79b487b3a7d395f)
127d8b95cSPierre Jolivet #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/
2ebda224cSfranklin5 
3ebda224cSfranklin5 typedef struct {
4ebda224cSfranklin5   Mat A;
5a48507d3SJose E. Roman   Mat D; /* local submatrix for diagonal part */
627d8b95cSPierre Jolivet   Vec w;
713b343c1SBarry Smith } Mat_NormalHermitian;
8ebda224cSfranklin5 
966976f2fSJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_NormalHermitian(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
10d71ae5a4SJacob Faibussowitsch {
1127d8b95cSPierre Jolivet   Mat_NormalHermitian *a;
1227d8b95cSPierre Jolivet   Mat                  B, *suba;
13ab704866SPierre Jolivet   IS                  *row;
14*b9c875b8SPierre Jolivet   PetscScalar          shift, scale;
15ab704866SPierre Jolivet   PetscInt             M;
16ab704866SPierre Jolivet 
17ab704866SPierre Jolivet   PetscFunctionBegin;
1827d8b95cSPierre Jolivet   PetscCheck(irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented");
19*b9c875b8SPierre Jolivet   PetscCall(MatShellGetScalingShifts(mat, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
2027d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(mat, &a));
2127d8b95cSPierre Jolivet   B = a->A;
2248a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
23ab704866SPierre Jolivet   PetscCall(MatGetSize(B, &M, NULL));
24ab704866SPierre Jolivet   PetscCall(PetscMalloc1(n, &row));
25ab704866SPierre Jolivet   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0]));
26ab704866SPierre Jolivet   PetscCall(ISSetIdentity(row[0]));
27ab704866SPierre Jolivet   for (M = 1; M < n; ++M) row[M] = row[0];
28ab704866SPierre Jolivet   PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba));
29ab704866SPierre Jolivet   for (M = 0; M < n; ++M) {
30ab704866SPierre Jolivet     PetscCall(MatCreateNormalHermitian(suba[M], *submat + M));
31*b9c875b8SPierre Jolivet     PetscCall(MatShift((*submat)[M], shift));
32*b9c875b8SPierre Jolivet     PetscCall(MatScale((*submat)[M], scale));
33ab704866SPierre Jolivet   }
34ab704866SPierre Jolivet   PetscCall(ISDestroy(&row[0]));
35ab704866SPierre Jolivet   PetscCall(PetscFree(row));
36ab704866SPierre Jolivet   PetscCall(MatDestroySubMatrices(n, &suba));
373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38ab704866SPierre Jolivet }
39ab704866SPierre Jolivet 
4066976f2fSJacob Faibussowitsch static PetscErrorCode MatPermute_NormalHermitian(Mat A, IS rowp, IS colp, Mat *B)
41d71ae5a4SJacob Faibussowitsch {
4227d8b95cSPierre Jolivet   Mat_NormalHermitian *a;
4327d8b95cSPierre Jolivet   Mat                  C, Aa;
44ab704866SPierre Jolivet   IS                   row;
45*b9c875b8SPierre Jolivet   PetscScalar          shift, scale;
46ab704866SPierre Jolivet 
47ab704866SPierre Jolivet   PetscFunctionBegin;
48ab704866SPierre Jolivet   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
49*b9c875b8SPierre Jolivet   PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
5027d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
5127d8b95cSPierre Jolivet   Aa = a->A;
52ab704866SPierre Jolivet   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row));
53ab704866SPierre Jolivet   PetscCall(ISSetIdentity(row));
54ab704866SPierre Jolivet   PetscCall(MatPermute(Aa, row, colp, &C));
55ab704866SPierre Jolivet   PetscCall(ISDestroy(&row));
56ab704866SPierre Jolivet   PetscCall(MatCreateNormalHermitian(C, B));
57ab704866SPierre Jolivet   PetscCall(MatDestroy(&C));
58*b9c875b8SPierre Jolivet   PetscCall(MatShift(*B, shift));
59*b9c875b8SPierre Jolivet   PetscCall(MatScale(*B, scale));
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61ab704866SPierre Jolivet }
62ab704866SPierre Jolivet 
6366976f2fSJacob Faibussowitsch static PetscErrorCode MatDuplicate_NormalHermitian(Mat A, MatDuplicateOption op, Mat *B)
64d71ae5a4SJacob Faibussowitsch {
6527d8b95cSPierre Jolivet   Mat_NormalHermitian *a;
66ab704866SPierre Jolivet   Mat                  C;
67ab704866SPierre Jolivet 
68ab704866SPierre Jolivet   PetscFunctionBegin;
6927d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
70ab704866SPierre Jolivet   PetscCall(MatDuplicate(a->A, op, &C));
71ab704866SPierre Jolivet   PetscCall(MatCreateNormalHermitian(C, B));
72ab704866SPierre Jolivet   PetscCall(MatDestroy(&C));
7327d8b95cSPierre Jolivet   if (op == MAT_COPY_VALUES) PetscCall(MatCopy(A, *B, SAME_NONZERO_PATTERN));
743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75ab704866SPierre Jolivet }
76ab704866SPierre Jolivet 
7766976f2fSJacob Faibussowitsch static PetscErrorCode MatCopy_NormalHermitian(Mat A, Mat B, MatStructure str)
78d71ae5a4SJacob Faibussowitsch {
7927d8b95cSPierre Jolivet   Mat_NormalHermitian *a, *b;
80ab704866SPierre Jolivet 
81ab704866SPierre Jolivet   PetscFunctionBegin;
8227d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
8327d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(B, &b));
84ab704866SPierre Jolivet   PetscCall(MatCopy(a->A, b->A, str));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86ab704866SPierre Jolivet }
87ab704866SPierre Jolivet 
8866976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_NormalHermitian(Mat N, Vec x, Vec y)
89d71ae5a4SJacob Faibussowitsch {
9027d8b95cSPierre Jolivet   Mat_NormalHermitian *Na;
91ebda224cSfranklin5 
92ebda224cSfranklin5   PetscFunctionBegin;
9327d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
9427d8b95cSPierre Jolivet   PetscCall(MatMult(Na->A, x, Na->w));
959566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(Na->A, Na->w, y));
963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97ebda224cSfranklin5 }
98ebda224cSfranklin5 
9966976f2fSJacob Faibussowitsch static PetscErrorCode MatDestroy_NormalHermitian(Mat N)
100d71ae5a4SJacob Faibussowitsch {
10127d8b95cSPierre Jolivet   Mat_NormalHermitian *Na;
102ebda224cSfranklin5 
103ebda224cSfranklin5   PetscFunctionBegin;
10427d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
1059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
106a48507d3SJose E. Roman   PetscCall(MatDestroy(&Na->D));
1079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->w));
10827d8b95cSPierre Jolivet   PetscCall(PetscFree(Na));
10914e4dea2SJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalHermitianGetMat_C", NULL));
11027d8b95cSPierre Jolivet #if !defined(PETSC_USE_COMPLEX)
11127d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMat_C", NULL));
11227d8b95cSPierre Jolivet #endif
113ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_seqaij_C", NULL));
114ab704866SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_mpiaij_C", NULL));
1152f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
1162f0a5c5aSJose E. Roman   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normalh_hypre_C", NULL));
1172f0a5c5aSJose E. Roman #endif
11827d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatShellSetContext_C", NULL));
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120ebda224cSfranklin5 }
121ebda224cSfranklin5 
122ebda224cSfranklin5 /*
123ebda224cSfranklin5       Slow, nonscalable version
124ebda224cSfranklin5 */
12566976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_NormalHermitian(Mat N, Vec v)
126d71ae5a4SJacob Faibussowitsch {
12727d8b95cSPierre Jolivet   Mat_NormalHermitian *Na;
12827d8b95cSPierre Jolivet   Mat                  A;
129ebda224cSfranklin5   PetscInt             i, j, rstart, rend, nnz;
130ebda224cSfranklin5   const PetscInt      *cols;
131ebda224cSfranklin5   PetscScalar         *diag, *work, *values;
132ebda224cSfranklin5   const PetscScalar   *mvalues;
133ebda224cSfranklin5 
134ebda224cSfranklin5   PetscFunctionBegin;
13527d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
13627d8b95cSPierre Jolivet   A = Na->A;
1379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work));
1389566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(work, A->cmap->N));
1399566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
140ebda224cSfranklin5   for (i = rstart; i < rend; i++) {
1419566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues));
142ad540459SPierre Jolivet     for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * PetscConj(mvalues[j]);
1439566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues));
144ebda224cSfranklin5   }
1451c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N)));
146ebda224cSfranklin5   rstart = N->cmap->rstart;
147ebda224cSfranklin5   rend   = N->cmap->rend;
1489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &values));
1499566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart));
1509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &values));
1519566063dSJacob Faibussowitsch   PetscCall(PetscFree2(diag, work));
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
153ebda224cSfranklin5 }
154ebda224cSfranklin5 
15566976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_NormalHermitian(Mat N, Mat *D)
156d71ae5a4SJacob Faibussowitsch {
15727d8b95cSPierre Jolivet   Mat_NormalHermitian *Na;
15827d8b95cSPierre Jolivet   Mat                  M, A;
159a48507d3SJose E. Roman 
160a48507d3SJose E. Roman   PetscFunctionBegin;
16127d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)N->data)->zrows && !((Mat_Shell *)N->data)->zcols, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME
16227d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)N->data)->axpy, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatAXPY() has been called on the input Mat");                                            // TODO FIXME
16327d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)N->data)->left && !((Mat_Shell *)N->data)->right, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME
16427d8b95cSPierre Jolivet   PetscCheck(!((Mat_Shell *)N->data)->dshift, PetscObjectComm((PetscObject)N), PETSC_ERR_SUP, "Cannot call MatGetDiagonalBlock() if MatDiagonalSet() has been called on the input Mat");                                   // TODO FIXME
16527d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(N, &Na));
16627d8b95cSPierre Jolivet   A = Na->A;
167a48507d3SJose E. Roman   PetscCall(MatGetDiagonalBlock(A, &M));
168a48507d3SJose E. Roman   PetscCall(MatCreateNormalHermitian(M, &Na->D));
169a48507d3SJose E. Roman   *D = Na->D;
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
171a48507d3SJose E. Roman }
172a48507d3SJose E. Roman 
17314e4dea2SJose E. Roman static PetscErrorCode MatNormalHermitianGetMat_NormalHermitian(Mat A, Mat *M)
174d71ae5a4SJacob Faibussowitsch {
17527d8b95cSPierre Jolivet   Mat_NormalHermitian *Aa;
17665f45395SPierre Jolivet 
17765f45395SPierre Jolivet   PetscFunctionBegin;
17827d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &Aa));
17965f45395SPierre Jolivet   *M = Aa->A;
1803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18165f45395SPierre Jolivet }
18265f45395SPierre Jolivet 
18365f45395SPierre Jolivet /*@
18411a5261eSBarry Smith   MatNormalHermitianGetMat - Gets the `Mat` object stored inside a `MATNORMALHERMITIAN`
18565f45395SPierre Jolivet 
1862ef1f0ffSBarry Smith   Logically Collective
18765f45395SPierre Jolivet 
18865f45395SPierre Jolivet   Input Parameter:
18911a5261eSBarry Smith . A - the `MATNORMALHERMITIAN` matrix
19065f45395SPierre Jolivet 
19165f45395SPierre Jolivet   Output Parameter:
19227d8b95cSPierre Jolivet . M - the matrix object stored inside `A`
19365f45395SPierre Jolivet 
19465f45395SPierre Jolivet   Level: intermediate
19565f45395SPierre Jolivet 
1961cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
19765f45395SPierre Jolivet @*/
198d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormalHermitianGetMat(Mat A, Mat *M)
199d71ae5a4SJacob Faibussowitsch {
20065f45395SPierre Jolivet   PetscFunctionBegin;
20165f45395SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
20265f45395SPierre Jolivet   PetscValidType(A, 1);
2034f572ea9SToby Isaac   PetscAssertPointer(M, 2);
20414e4dea2SJose E. Roman   PetscUseMethod(A, "MatNormalHermitianGetMat_C", (Mat, Mat *), (A, M));
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20665f45395SPierre Jolivet }
20765f45395SPierre Jolivet 
20866976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_NormalHermitian_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
209d71ae5a4SJacob Faibussowitsch {
21027d8b95cSPierre Jolivet   Mat_NormalHermitian *Aa;
211ab704866SPierre Jolivet   Mat                  B, conjugate;
212ab704866SPierre Jolivet   PetscInt             m, n, M, N;
213ab704866SPierre Jolivet 
214ab704866SPierre Jolivet   PetscFunctionBegin;
21527d8b95cSPierre Jolivet   PetscCall(MatShellGetContext(A, &Aa));
216ab704866SPierre Jolivet   PetscCall(MatGetSize(A, &M, &N));
217ab704866SPierre Jolivet   PetscCall(MatGetLocalSize(A, &m, &n));
218ab704866SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
219ab704866SPierre Jolivet     B = *newmat;
220ab704866SPierre Jolivet     PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B));
221ab704866SPierre Jolivet   } else {
222ab704866SPierre Jolivet     PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B));
223ab704866SPierre Jolivet     PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
224ab704866SPierre Jolivet     PetscCall(MatProductSetFromOptions(B));
225ab704866SPierre Jolivet     PetscCall(MatProductSymbolic(B));
226ab704866SPierre Jolivet     PetscCall(MatSetOption(B, !PetscDefined(USE_COMPLEX) ? MAT_SYMMETRIC : MAT_HERMITIAN, PETSC_TRUE));
227ab704866SPierre Jolivet   }
228ab704866SPierre Jolivet   if (PetscDefined(USE_COMPLEX)) {
229ab704866SPierre Jolivet     PetscCall(MatDuplicate(Aa->A, MAT_COPY_VALUES, &conjugate));
230ab704866SPierre Jolivet     PetscCall(MatConjugate(conjugate));
231ab704866SPierre Jolivet     PetscCall(MatProductReplaceMats(conjugate, Aa->A, NULL, B));
232ab704866SPierre Jolivet   }
233ab704866SPierre Jolivet   PetscCall(MatProductNumeric(B));
234ab704866SPierre Jolivet   if (PetscDefined(USE_COMPLEX)) PetscCall(MatDestroy(&conjugate));
235ab704866SPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
236ab704866SPierre Jolivet     PetscCall(MatHeaderReplace(A, &B));
237ab704866SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
238ab704866SPierre Jolivet   PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat));
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240ab704866SPierre Jolivet }
241ab704866SPierre Jolivet 
2422f0a5c5aSJose E. Roman #if defined(PETSC_HAVE_HYPRE)
24366976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_NormalHermitian_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
2442f0a5c5aSJose E. Roman {
2452f0a5c5aSJose E. Roman   PetscFunctionBegin;
2462f0a5c5aSJose E. Roman   if (reuse == MAT_INITIAL_MATRIX) {
2472f0a5c5aSJose E. Roman     PetscCall(MatConvert(A, MATAIJ, reuse, B));
2482f0a5c5aSJose E. Roman     PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B));
2492f0a5c5aSJose E. Roman   } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2512f0a5c5aSJose E. Roman }
2522f0a5c5aSJose E. Roman #endif
2532f0a5c5aSJose E. Roman 
2542ef1f0ffSBarry Smith /*MC
2552ef1f0ffSBarry Smith   MATNORMALHERMITIAN - a matrix that behaves like (A*)'*A for `MatMult()` while only containing A
2562ef1f0ffSBarry Smith 
2572ef1f0ffSBarry Smith   Level: intermediate
2582ef1f0ffSBarry Smith 
25927d8b95cSPierre Jolivet   Developer Notes:
26027d8b95cSPierre Jolivet   This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
26127d8b95cSPierre Jolivet 
26227d8b95cSPierre Jolivet   Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage
26327d8b95cSPierre Jolivet 
2641cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateNormalHermitian()`, `MatMult()`, `MatNormalHermitianGetMat()`, `MATNORMAL`, `MatCreateNormal()`
2652ef1f0ffSBarry Smith M*/
2662ef1f0ffSBarry Smith 
267ebda224cSfranklin5 /*@
26811a5261eSBarry Smith   MatCreateNormalHermitian - Creates a new matrix object `MATNORMALHERMITIAN` that behaves like (A*)'*A.
269ebda224cSfranklin5 
270c3339decSBarry Smith   Collective
271ebda224cSfranklin5 
272ebda224cSfranklin5   Input Parameter:
273ebda224cSfranklin5 . A - the (possibly rectangular complex) matrix
274ebda224cSfranklin5 
275ebda224cSfranklin5   Output Parameter:
276ebda224cSfranklin5 . N - the matrix that represents (A*)'*A
277ebda224cSfranklin5 
278ebda224cSfranklin5   Level: intermediate
279ebda224cSfranklin5 
28011a5261eSBarry Smith   Note:
28195452b02SPatrick Sanan   The product (A*)'*A is NOT actually formed! Rather the new matrix
28211a5261eSBarry Smith   object performs the matrix-vector product, `MatMult()`, by first multiplying by
283ebda224cSfranklin5   A and then (A*)'
28411a5261eSBarry Smith 
2851cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatNormalHermitianGetMat()`
286ebda224cSfranklin5 @*/
287d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNormalHermitian(Mat A, Mat *N)
288d71ae5a4SJacob Faibussowitsch {
28913b343c1SBarry Smith   Mat_NormalHermitian *Na;
2905cb049f5SJose E. Roman   VecType              vtype;
291ebda224cSfranklin5 
292ebda224cSfranklin5   PetscFunctionBegin;
2939566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N));
2949566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap));
2959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap));
29627d8b95cSPierre Jolivet   PetscCall(MatSetType(*N, MATSHELL));
2974dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&Na));
29827d8b95cSPierre Jolivet   PetscCall(MatShellSetContext(*N, Na));
2999566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
300ebda224cSfranklin5   Na->A = A;
3019566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, NULL, &Na->w));
302ebda224cSfranklin5 
30327d8b95cSPierre Jolivet   PetscCall(MatSetBlockSizes(*N, PetscAbs(A->cmap->bs), PetscAbs(A->rmap->bs)));
30427d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_DESTROY, (void (*)(void))MatDestroy_NormalHermitian));
30527d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_MULT, (void (*)(void))MatMult_NormalHermitian));
30627d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_MULT_HERMITIAN_TRANSPOSE, (void (*)(void))MatMult_NormalHermitian));
30727d8b95cSPierre Jolivet #if !defined(PETSC_USE_COMPLEX)
30827d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_NormalHermitian));
30927d8b95cSPierre Jolivet #endif
31027d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_NormalHermitian));
31127d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NormalHermitian));
31227d8b95cSPierre Jolivet   PetscCall(MatShellSetOperation(*N, MATOP_COPY, (void (*)(void))MatCopy_NormalHermitian));
3132e5cc420SJose E. Roman   (*N)->ops->getdiagonalblock  = MatGetDiagonalBlock_NormalHermitian;
314ab704866SPierre Jolivet   (*N)->ops->createsubmatrices = MatCreateSubMatrices_NormalHermitian;
315ab704866SPierre Jolivet   (*N)->ops->permute           = MatPermute_NormalHermitian;
31612ab1b64SPierre Jolivet 
31727d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalHermitianGetMat_C", MatNormalHermitianGetMat_NormalHermitian));
31827d8b95cSPierre Jolivet #if !defined(PETSC_USE_COMPLEX)
31927d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalGetMat_C", MatNormalHermitianGetMat_NormalHermitian));
3202f0a5c5aSJose E. Roman #endif
32127d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normalh_seqaij_C", MatConvert_NormalHermitian_AIJ));
32227d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normalh_mpiaij_C", MatConvert_NormalHermitian_AIJ));
32327d8b95cSPierre Jolivet #if defined(PETSC_HAVE_HYPRE)
32427d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normalh_hypre_C", MatConvert_NormalHermitian_HYPRE));
32527d8b95cSPierre Jolivet #endif
32627d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContext_C", MatShellSetContext_Immutable));
32727d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
32827d8b95cSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
32927d8b95cSPierre Jolivet   PetscCall(MatSetOption(*N, !PetscDefined(USE_COMPLEX) ? MAT_SYMMETRIC : MAT_HERMITIAN, PETSC_TRUE));
3309566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A, &vtype));
3319566063dSJacob Faibussowitsch   PetscCall(MatSetVecType(*N, vtype));
3325cb049f5SJose E. Roman #if defined(PETSC_HAVE_DEVICE)
3339566063dSJacob Faibussowitsch   PetscCall(MatBindToCPU(*N, A->boundtocpu));
3345cb049f5SJose E. Roman #endif
33527d8b95cSPierre Jolivet   PetscCall(MatSetUp(*N));
33627d8b95cSPierre Jolivet   PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMALHERMITIAN));
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
338ebda224cSfranklin5 }
339