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