xref: /petsc/src/mat/impls/normal/normm.c (revision 4ad8454beace47809662cdae21ee081016eaa39a)
1 #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/
2 
3 typedef struct {
4   Mat A;
5   Mat D; /* local submatrix for diagonal part */
6   Vec w;
7 } Mat_Normal;
8 
9 static PetscErrorCode MatIncreaseOverlap_Normal(Mat A, PetscInt is_max, IS is[], PetscInt ov)
10 {
11   Mat_Normal *a;
12   Mat         pattern;
13 
14   PetscFunctionBegin;
15   PetscCheck(ov >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified");
16   PetscCall(MatShellGetContext(A, &a));
17   PetscCall(MatProductCreate(a->A, a->A, NULL, &pattern));
18   PetscCall(MatProductSetType(pattern, MATPRODUCT_AtB));
19   PetscCall(MatProductSetFromOptions(pattern));
20   PetscCall(MatProductSymbolic(pattern));
21   PetscCall(MatIncreaseOverlap(pattern, is_max, is, ov));
22   PetscCall(MatDestroy(&pattern));
23   PetscFunctionReturn(PETSC_SUCCESS);
24 }
25 
26 static PetscErrorCode MatCreateSubMatrices_Normal(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
27 {
28   Mat_Normal *a;
29   Mat         B, *suba;
30   IS         *row;
31   PetscInt    M;
32 
33   PetscFunctionBegin;
34   PetscCheck(!((Mat_Shell *)mat->data)->zrows && !((Mat_Shell *)mat->data)->zcols, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatZeroRows()/MatZeroRowsColumns() after the SubMatrices creation
35   PetscCheck(!((Mat_Shell *)mat->data)->axpy, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatAXPY() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatAXPY() after the SubMatrices creation
36   PetscCheck(!((Mat_Shell *)mat->data)->left && !((Mat_Shell *)mat->data)->right, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalScale() after the SubMatrices creation with a SubVector
37   PetscCheck(!((Mat_Shell *)mat->data)->dshift, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalSet() after the SubMatrices creation with a SubVector
38   PetscCheck(irow == icol, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Not implemented");
39   PetscCall(MatShellGetContext(mat, &a));
40   B = a->A;
41   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
42   PetscCall(MatGetSize(B, &M, NULL));
43   PetscCall(PetscMalloc1(n, &row));
44   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &row[0]));
45   PetscCall(ISSetIdentity(row[0]));
46   for (M = 1; M < n; ++M) row[M] = row[0];
47   PetscCall(MatCreateSubMatrices(B, n, row, icol, MAT_INITIAL_MATRIX, &suba));
48   for (M = 0; M < n; ++M) {
49     PetscCall(MatCreateNormal(suba[M], *submat + M));
50     ((Mat_Shell *)(*submat)[M]->data)->vscale = ((Mat_Shell *)mat->data)->vscale;
51     ((Mat_Shell *)(*submat)[M]->data)->vshift = ((Mat_Shell *)mat->data)->vshift;
52   }
53   PetscCall(ISDestroy(&row[0]));
54   PetscCall(PetscFree(row));
55   PetscCall(MatDestroySubMatrices(n, &suba));
56   PetscFunctionReturn(PETSC_SUCCESS);
57 }
58 
59 static PetscErrorCode MatPermute_Normal(Mat A, IS rowp, IS colp, Mat *B)
60 {
61   Mat_Normal *a;
62   Mat         C, Aa;
63   IS          row;
64 
65   PetscFunctionBegin;
66   PetscCheck(!((Mat_Shell *)A->data)->zrows && !((Mat_Shell *)A->data)->zcols, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatPermute() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatZeroRows()/MatZeroRowsColumns() after the permutation with a permuted zrows and zcols
67   PetscCheck(!((Mat_Shell *)A->data)->axpy, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatPermute() if MatAXPY() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatAXPY() after the permutation with a permuted axpy
68   PetscCheck(!((Mat_Shell *)A->data)->left && !((Mat_Shell *)A->data)->right, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatPermute() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalScale() after the permutation with a permuted left and right
69   PetscCheck(!((Mat_Shell *)A->data)->dshift, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatPermute() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME: lift this limitation by calling MatDiagonalSet() after the permutation with a permuted dshift
70   PetscCheck(rowp == colp, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Row permutation and column permutation must be the same");
71   PetscCall(MatShellGetContext(A, &a));
72   Aa = a->A;
73   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa), Aa->rmap->n, Aa->rmap->rstart, 1, &row));
74   PetscCall(ISSetIdentity(row));
75   PetscCall(MatPermute(Aa, row, colp, &C));
76   PetscCall(ISDestroy(&row));
77   PetscCall(MatCreateNormal(C, B));
78   PetscCall(MatDestroy(&C));
79   ((Mat_Shell *)(*B)->data)->vscale = ((Mat_Shell *)A->data)->vscale;
80   ((Mat_Shell *)(*B)->data)->vshift = ((Mat_Shell *)A->data)->vshift;
81   PetscFunctionReturn(PETSC_SUCCESS);
82 }
83 
84 static PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B)
85 {
86   Mat_Normal *a;
87   Mat         C;
88 
89   PetscFunctionBegin;
90   PetscCall(MatShellGetContext(A, &a));
91   PetscCall(MatDuplicate(a->A, op, &C));
92   PetscCall(MatCreateNormal(C, B));
93   PetscCall(MatDestroy(&C));
94   if (op == MAT_COPY_VALUES) PetscCall(MatCopy(A, *B, SAME_NONZERO_PATTERN));
95   PetscFunctionReturn(PETSC_SUCCESS);
96 }
97 
98 static PetscErrorCode MatCopy_Normal(Mat A, Mat B, MatStructure str)
99 {
100   Mat_Normal *a, *b;
101 
102   PetscFunctionBegin;
103   PetscCall(MatShellGetContext(A, &a));
104   PetscCall(MatShellGetContext(B, &b));
105   PetscCall(MatCopy(a->A, b->A, str));
106   PetscFunctionReturn(PETSC_SUCCESS);
107 }
108 
109 static PetscErrorCode MatMult_Normal(Mat N, Vec x, Vec y)
110 {
111   Mat_Normal *Na;
112 
113   PetscFunctionBegin;
114   PetscCall(MatShellGetContext(N, &Na));
115   PetscCall(MatMult(Na->A, x, Na->w));
116   PetscCall(MatMultTranspose(Na->A, Na->w, y));
117   PetscFunctionReturn(PETSC_SUCCESS);
118 }
119 
120 static PetscErrorCode MatDestroy_Normal(Mat N)
121 {
122   Mat_Normal *Na;
123 
124   PetscFunctionBegin;
125   PetscCall(MatShellGetContext(N, &Na));
126   PetscCall(MatDestroy(&Na->A));
127   PetscCall(MatDestroy(&Na->D));
128   PetscCall(VecDestroy(&Na->w));
129   PetscCall(PetscFree(Na));
130   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatNormalGetMat_C", NULL));
131   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_seqaij_C", NULL));
132   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_mpiaij_C", NULL));
133 #if defined(PETSC_HAVE_HYPRE)
134   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatConvert_normal_hypre_C", NULL));
135 #endif
136   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_seqdense_C", NULL));
137   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_mpidense_C", NULL));
138   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_normal_dense_C", NULL));
139   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatShellSetContext_C", NULL));
140   PetscFunctionReturn(PETSC_SUCCESS);
141 }
142 
143 /*
144       Slow, nonscalable version
145 */
146 static PetscErrorCode MatGetDiagonal_Normal(Mat N, Vec v)
147 {
148   Mat_Normal        *Na;
149   Mat                A;
150   PetscInt           i, j, rstart, rend, nnz;
151   const PetscInt    *cols;
152   PetscScalar       *diag, *work, *values;
153   const PetscScalar *mvalues;
154 
155   PetscFunctionBegin;
156   PetscCall(MatShellGetContext(N, &Na));
157   A = Na->A;
158   PetscCall(PetscMalloc2(A->cmap->N, &diag, A->cmap->N, &work));
159   PetscCall(PetscArrayzero(work, A->cmap->N));
160   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
161   for (i = rstart; i < rend; i++) {
162     PetscCall(MatGetRow(A, i, &nnz, &cols, &mvalues));
163     for (j = 0; j < nnz; j++) work[cols[j]] += mvalues[j] * mvalues[j];
164     PetscCall(MatRestoreRow(A, i, &nnz, &cols, &mvalues));
165   }
166   PetscCall(MPIU_Allreduce(work, diag, A->cmap->N, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)N)));
167   rstart = N->cmap->rstart;
168   rend   = N->cmap->rend;
169   PetscCall(VecGetArray(v, &values));
170   PetscCall(PetscArraycpy(values, diag + rstart, rend - rstart));
171   PetscCall(VecRestoreArray(v, &values));
172   PetscCall(PetscFree2(diag, work));
173   PetscFunctionReturn(PETSC_SUCCESS);
174 }
175 
176 static PetscErrorCode MatGetDiagonalBlock_Normal(Mat N, Mat *D)
177 {
178   Mat_Normal *Na;
179   Mat         M, A;
180 
181   PetscFunctionBegin;
182   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
183   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
184   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
185   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
186   PetscCall(MatShellGetContext(N, &Na));
187   A = Na->A;
188   PetscCall(MatGetDiagonalBlock(A, &M));
189   PetscCall(MatCreateNormal(M, &Na->D));
190   *D = Na->D;
191   PetscFunctionReturn(PETSC_SUCCESS);
192 }
193 
194 static PetscErrorCode MatNormalGetMat_Normal(Mat A, Mat *M)
195 {
196   Mat_Normal *Aa;
197 
198   PetscFunctionBegin;
199   PetscCall(MatShellGetContext(A, &Aa));
200   *M = Aa->A;
201   PetscFunctionReturn(PETSC_SUCCESS);
202 }
203 
204 /*@
205   MatNormalGetMat - Gets the `Mat` object stored inside a `MATNORMAL`
206 
207   Logically Collective
208 
209   Input Parameter:
210 . A - the `MATNORMAL` matrix
211 
212   Output Parameter:
213 . M - the matrix object stored inside `A`
214 
215   Level: intermediate
216 
217 .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatCreateNormal()`
218 @*/
219 PetscErrorCode MatNormalGetMat(Mat A, Mat *M)
220 {
221   PetscFunctionBegin;
222   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
223   PetscValidType(A, 1);
224   PetscAssertPointer(M, 2);
225   PetscUseMethod(A, "MatNormalGetMat_C", (Mat, Mat *), (A, M));
226   PetscFunctionReturn(PETSC_SUCCESS);
227 }
228 
229 static PetscErrorCode MatConvert_Normal_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
230 {
231   Mat_Normal *Aa;
232   Mat         B;
233   PetscInt    m, n, M, N;
234 
235   PetscFunctionBegin;
236   PetscCall(MatShellGetContext(A, &Aa));
237   PetscCall(MatGetSize(A, &M, &N));
238   PetscCall(MatGetLocalSize(A, &m, &n));
239   if (reuse == MAT_REUSE_MATRIX) {
240     B = *newmat;
241     PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B));
242   } else {
243     PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B));
244     PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
245     PetscCall(MatProductSetFromOptions(B));
246     PetscCall(MatProductSymbolic(B));
247     PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE));
248   }
249   PetscCall(MatProductNumeric(B));
250   if (reuse == MAT_INPLACE_MATRIX) {
251     PetscCall(MatHeaderReplace(A, &B));
252   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
253   PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat));
254   PetscFunctionReturn(PETSC_SUCCESS);
255 }
256 
257 #if defined(PETSC_HAVE_HYPRE)
258 static PetscErrorCode MatConvert_Normal_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
259 {
260   PetscFunctionBegin;
261   if (reuse == MAT_INITIAL_MATRIX) {
262     PetscCall(MatConvert(A, MATAIJ, reuse, B));
263     PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B));
264   } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */
265   PetscFunctionReturn(PETSC_SUCCESS);
266 }
267 #endif
268 
269 typedef struct {
270   Mat work[2];
271 } Normal_Dense;
272 
273 static PetscErrorCode MatProductNumeric_Normal_Dense(Mat C)
274 {
275   Mat           A, B;
276   Normal_Dense *contents;
277   Mat_Normal   *a;
278   PetscScalar  *array;
279 
280   PetscFunctionBegin;
281   MatCheckProduct(C, 1);
282   A = C->product->A;
283   B = C->product->B;
284   PetscCall(MatShellGetContext(A, &a));
285   contents = (Normal_Dense *)C->product->data;
286   PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
287   if (((Mat_Shell *)A->data)->right) {
288     PetscCall(MatCopy(B, C, SAME_NONZERO_PATTERN));
289     PetscCall(MatDiagonalScale(C, ((Mat_Shell *)A->data)->right, NULL));
290   }
291   PetscCall(MatProductNumeric(contents->work[0]));
292   PetscCall(MatDenseGetArrayWrite(C, &array));
293   PetscCall(MatDensePlaceArray(contents->work[1], array));
294   PetscCall(MatProductNumeric(contents->work[1]));
295   PetscCall(MatDenseRestoreArrayWrite(C, &array));
296   PetscCall(MatDenseResetArray(contents->work[1]));
297   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
298   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
299   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
300   PetscCall(MatScale(C, ((Mat_Shell *)A->data)->vscale));
301   PetscFunctionReturn(PETSC_SUCCESS);
302 }
303 
304 static PetscErrorCode MatNormal_DenseDestroy(void *ctx)
305 {
306   Normal_Dense *contents = (Normal_Dense *)ctx;
307 
308   PetscFunctionBegin;
309   PetscCall(MatDestroy(contents->work));
310   PetscCall(MatDestroy(contents->work + 1));
311   PetscCall(PetscFree(contents));
312   PetscFunctionReturn(PETSC_SUCCESS);
313 }
314 
315 static PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C)
316 {
317   Mat           A, B;
318   Normal_Dense *contents = NULL;
319   Mat_Normal   *a;
320   PetscScalar  *array;
321   PetscInt      n, N, m, M;
322 
323   PetscFunctionBegin;
324   MatCheckProduct(C, 1);
325   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
326   A = C->product->A;
327   B = C->product->B;
328   PetscCall(MatShellGetContext(A, &a));
329   PetscCheck(!((Mat_Shell *)A->data)->zrows && !((Mat_Shell *)A->data)->zcols, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatProductSymbolic() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME
330   PetscCheck(!((Mat_Shell *)A->data)->axpy, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatProductSymbolic() if MatAXPY() has been called on the input Mat");          // TODO FIXME
331   PetscCheck(!((Mat_Shell *)A->data)->left, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatProductSymbolic() if MatDiagonalScale() has been called on the input Mat"); // TODO FIXME
332   PetscCheck(!((Mat_Shell *)A->data)->dshift, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatProductSymbolic() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME
333   PetscCall(MatGetLocalSize(C, &m, &n));
334   PetscCall(MatGetSize(C, &M, &N));
335   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
336     PetscCall(MatGetLocalSize(B, NULL, &n));
337     PetscCall(MatGetSize(B, NULL, &N));
338     PetscCall(MatGetLocalSize(A, &m, NULL));
339     PetscCall(MatGetSize(A, &M, NULL));
340     PetscCall(MatSetSizes(C, m, n, M, N));
341   }
342   PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
343   PetscCall(MatSetUp(C));
344   PetscCall(PetscNew(&contents));
345   C->product->data    = contents;
346   C->product->destroy = MatNormal_DenseDestroy;
347   if (((Mat_Shell *)A->data)->right) {
348     PetscCall(MatProductCreate(a->A, C, NULL, contents->work));
349   } else {
350     PetscCall(MatProductCreate(a->A, B, NULL, contents->work));
351   }
352   PetscCall(MatProductSetType(contents->work[0], MATPRODUCT_AB));
353   PetscCall(MatProductSetFromOptions(contents->work[0]));
354   PetscCall(MatProductSymbolic(contents->work[0]));
355   PetscCall(MatProductCreate(a->A, contents->work[0], NULL, contents->work + 1));
356   PetscCall(MatProductSetType(contents->work[1], MATPRODUCT_AtB));
357   PetscCall(MatProductSetFromOptions(contents->work[1]));
358   PetscCall(MatProductSymbolic(contents->work[1]));
359   PetscCall(MatDenseGetArrayWrite(C, &array));
360   PetscCall(MatSeqDenseSetPreallocation(contents->work[1], array));
361   PetscCall(MatMPIDenseSetPreallocation(contents->work[1], array));
362   PetscCall(MatDenseRestoreArrayWrite(C, &array));
363   C->ops->productnumeric = MatProductNumeric_Normal_Dense;
364   PetscFunctionReturn(PETSC_SUCCESS);
365 }
366 
367 static PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C)
368 {
369   PetscFunctionBegin;
370   C->ops->productsymbolic = MatProductSymbolic_Normal_Dense;
371   PetscFunctionReturn(PETSC_SUCCESS);
372 }
373 
374 static PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C)
375 {
376   Mat_Product *product = C->product;
377 
378   PetscFunctionBegin;
379   if (product->type == MATPRODUCT_AB) PetscCall(MatProductSetFromOptions_Normal_Dense_AB(C));
380   PetscFunctionReturn(PETSC_SUCCESS);
381 }
382 
383 /*MC
384   MATNORMAL - a matrix that behaves like A'*A for `MatMult()` while only containing A
385 
386   Level: intermediate
387 
388   Developer Notes:
389   This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
390 
391   Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage
392 
393 .seealso: [](ch_matrices), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
394 M*/
395 
396 /*@
397   MatCreateNormal - Creates a new `MATNORMAL` matrix object that behaves like A'*A.
398 
399   Collective
400 
401   Input Parameter:
402 . A - the (possibly rectangular) matrix
403 
404   Output Parameter:
405 . N - the matrix that represents A'*A
406 
407   Level: intermediate
408 
409   Notes:
410   The product A'*A is NOT actually formed! Rather the new matrix
411   object performs the matrix-vector product, `MatMult()`, by first multiplying by
412   A and then A'
413 
414 .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()`
415 @*/
416 PetscErrorCode MatCreateNormal(Mat A, Mat *N)
417 {
418   Mat_Normal *Na;
419   VecType     vtype;
420 
421   PetscFunctionBegin;
422   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N));
423   PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap));
424   PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap));
425   PetscCall(MatSetType(*N, MATSHELL));
426   PetscCall(PetscNew(&Na));
427   PetscCall(MatShellSetContext(*N, Na));
428   PetscCall(PetscObjectReference((PetscObject)A));
429   Na->A = A;
430   PetscCall(MatCreateVecs(A, NULL, &Na->w));
431 
432   PetscCall(MatSetBlockSizes(*N, PetscAbs(A->cmap->bs), PetscAbs(A->rmap->bs)));
433   PetscCall(MatShellSetOperation(*N, MATOP_DESTROY, (void (*)(void))MatDestroy_Normal));
434   PetscCall(MatShellSetOperation(*N, MATOP_MULT, (void (*)(void))MatMult_Normal));
435   PetscCall(MatShellSetOperation(*N, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_Normal));
436   PetscCall(MatShellSetOperation(*N, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_Normal));
437   PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Normal));
438   PetscCall(MatShellSetOperation(*N, MATOP_COPY, (void (*)(void))MatCopy_Normal));
439   (*N)->ops->getdiagonalblock  = MatGetDiagonalBlock_Normal;
440   (*N)->ops->increaseoverlap   = MatIncreaseOverlap_Normal;
441   (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal;
442   (*N)->ops->permute           = MatPermute_Normal;
443 
444   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalGetMat_C", MatNormalGetMat_Normal));
445   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_seqaij_C", MatConvert_Normal_AIJ));
446   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_mpiaij_C", MatConvert_Normal_AIJ));
447 #if defined(PETSC_HAVE_HYPRE)
448   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_hypre_C", MatConvert_Normal_HYPRE));
449 #endif
450   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_seqdense_C", MatProductSetFromOptions_Normal_Dense));
451   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_mpidense_C", MatProductSetFromOptions_Normal_Dense));
452   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_dense_C", MatProductSetFromOptions_Normal_Dense));
453   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContext_C", MatShellSetContext_Immutable));
454   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
455   PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
456   PetscCall(MatSetOption(*N, MAT_SYMMETRIC, PETSC_TRUE));
457   PetscCall(MatGetVecType(A, &vtype));
458   PetscCall(MatSetVecType(*N, vtype));
459 #if defined(PETSC_HAVE_DEVICE)
460   PetscCall(MatBindToCPU(*N, A->boundtocpu));
461 #endif
462   PetscCall(MatSetUp(*N));
463   PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMAL));
464   PetscFunctionReturn(PETSC_SUCCESS);
465 }
466