xref: /petsc/src/mat/impls/normal/normm.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 
2 #include <petsc/private/matimpl.h>          /*I "petscmat.h" I*/
3 
4 typedef struct {
5   Mat         A;
6   Vec         w,left,right,leftwork,rightwork;
7   PetscScalar scale;
8 } Mat_Normal;
9 
10 PetscErrorCode MatScale_Normal(Mat inA,PetscScalar scale)
11 {
12   Mat_Normal *a = (Mat_Normal*)inA->data;
13 
14   PetscFunctionBegin;
15   a->scale *= scale;
16   PetscFunctionReturn(0);
17 }
18 
19 PetscErrorCode MatDiagonalScale_Normal(Mat inA,Vec left,Vec right)
20 {
21   Mat_Normal     *a = (Mat_Normal*)inA->data;
22 
23   PetscFunctionBegin;
24   if (left) {
25     if (!a->left) {
26       CHKERRQ(VecDuplicate(left,&a->left));
27       CHKERRQ(VecCopy(left,a->left));
28     } else {
29       CHKERRQ(VecPointwiseMult(a->left,left,a->left));
30     }
31   }
32   if (right) {
33     if (!a->right) {
34       CHKERRQ(VecDuplicate(right,&a->right));
35       CHKERRQ(VecCopy(right,a->right));
36     } else {
37       CHKERRQ(VecPointwiseMult(a->right,right,a->right));
38     }
39   }
40   PetscFunctionReturn(0);
41 }
42 
43 PetscErrorCode MatIncreaseOverlap_Normal(Mat A,PetscInt is_max,IS is[],PetscInt ov)
44 {
45   Mat_Normal     *a = (Mat_Normal*)A->data;
46   Mat            pattern;
47 
48   PetscFunctionBegin;
49   PetscCheckFalse(ov < 0,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
50   CHKERRQ(MatProductCreate(a->A,a->A,NULL,&pattern));
51   CHKERRQ(MatProductSetType(pattern,MATPRODUCT_AtB));
52   CHKERRQ(MatProductSetFromOptions(pattern));
53   CHKERRQ(MatProductSymbolic(pattern));
54   CHKERRQ(MatIncreaseOverlap(pattern,is_max,is,ov));
55   CHKERRQ(MatDestroy(&pattern));
56   PetscFunctionReturn(0);
57 }
58 
59 PetscErrorCode MatCreateSubMatrices_Normal(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
60 {
61   Mat_Normal     *a = (Mat_Normal*)mat->data;
62   Mat            B = a->A, *suba;
63   IS             *row;
64   PetscInt       M;
65 
66   PetscFunctionBegin;
67   PetscCheckFalse(a->left || a->right || irow != icol,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not implemented");
68   if (scall != MAT_REUSE_MATRIX) {
69     CHKERRQ(PetscCalloc1(n,submat));
70   }
71   CHKERRQ(MatGetSize(B,&M,NULL));
72   CHKERRQ(PetscMalloc1(n,&row));
73   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,M,0,1,&row[0]));
74   CHKERRQ(ISSetIdentity(row[0]));
75   for (M = 1; M < n; ++M) row[M] = row[0];
76   CHKERRQ(MatCreateSubMatrices(B,n,row,icol,MAT_INITIAL_MATRIX,&suba));
77   for (M = 0; M < n; ++M) {
78     CHKERRQ(MatCreateNormal(suba[M],*submat+M));
79     ((Mat_Normal*)(*submat)[M]->data)->scale = a->scale;
80   }
81   CHKERRQ(ISDestroy(&row[0]));
82   CHKERRQ(PetscFree(row));
83   CHKERRQ(MatDestroySubMatrices(n,&suba));
84   PetscFunctionReturn(0);
85 }
86 
87 PetscErrorCode MatPermute_Normal(Mat A,IS rowp,IS colp,Mat *B)
88 {
89   Mat_Normal     *a = (Mat_Normal*)A->data;
90   Mat            C,Aa = a->A;
91   IS             row;
92 
93   PetscFunctionBegin;
94   PetscCheckFalse(rowp != colp,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Row permutation and column permutation must be the same");
95   CHKERRQ(ISCreateStride(PetscObjectComm((PetscObject)Aa),Aa->rmap->n,Aa->rmap->rstart,1,&row));
96   CHKERRQ(ISSetIdentity(row));
97   CHKERRQ(MatPermute(Aa,row,colp,&C));
98   CHKERRQ(ISDestroy(&row));
99   CHKERRQ(MatCreateNormal(C,B));
100   CHKERRQ(MatDestroy(&C));
101   PetscFunctionReturn(0);
102 }
103 
104 PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B)
105 {
106   Mat_Normal     *a = (Mat_Normal*)A->data;
107   Mat            C;
108 
109   PetscFunctionBegin;
110   PetscCheckFalse(a->left || a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented");
111   CHKERRQ(MatDuplicate(a->A,op,&C));
112   CHKERRQ(MatCreateNormal(C,B));
113   CHKERRQ(MatDestroy(&C));
114   if (op == MAT_COPY_VALUES) ((Mat_Normal*)(*B)->data)->scale = a->scale;
115   PetscFunctionReturn(0);
116 }
117 
118 PetscErrorCode MatCopy_Normal(Mat A,Mat B,MatStructure str)
119 {
120   Mat_Normal     *a = (Mat_Normal*)A->data,*b = (Mat_Normal*)B->data;
121 
122   PetscFunctionBegin;
123   PetscCheckFalse(a->left || a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented");
124   CHKERRQ(MatCopy(a->A,b->A,str));
125   b->scale = a->scale;
126   CHKERRQ(VecDestroy(&b->left));
127   CHKERRQ(VecDestroy(&b->right));
128   CHKERRQ(VecDestroy(&b->leftwork));
129   CHKERRQ(VecDestroy(&b->rightwork));
130   PetscFunctionReturn(0);
131 }
132 
133 PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y)
134 {
135   Mat_Normal     *Na = (Mat_Normal*)N->data;
136   Vec            in;
137 
138   PetscFunctionBegin;
139   in = x;
140   if (Na->right) {
141     if (!Na->rightwork) {
142       CHKERRQ(VecDuplicate(Na->right,&Na->rightwork));
143     }
144     CHKERRQ(VecPointwiseMult(Na->rightwork,Na->right,in));
145     in   = Na->rightwork;
146   }
147   CHKERRQ(MatMult(Na->A,in,Na->w));
148   CHKERRQ(MatMultTranspose(Na->A,Na->w,y));
149   if (Na->left) {
150     CHKERRQ(VecPointwiseMult(y,Na->left,y));
151   }
152   CHKERRQ(VecScale(y,Na->scale));
153   PetscFunctionReturn(0);
154 }
155 
156 PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
157 {
158   Mat_Normal     *Na = (Mat_Normal*)N->data;
159   Vec            in;
160 
161   PetscFunctionBegin;
162   in = v1;
163   if (Na->right) {
164     if (!Na->rightwork) {
165       CHKERRQ(VecDuplicate(Na->right,&Na->rightwork));
166     }
167     CHKERRQ(VecPointwiseMult(Na->rightwork,Na->right,in));
168     in   = Na->rightwork;
169   }
170   CHKERRQ(MatMult(Na->A,in,Na->w));
171   CHKERRQ(VecScale(Na->w,Na->scale));
172   if (Na->left) {
173     CHKERRQ(MatMultTranspose(Na->A,Na->w,v3));
174     CHKERRQ(VecPointwiseMult(v3,Na->left,v3));
175     CHKERRQ(VecAXPY(v3,1.0,v2));
176   } else {
177     CHKERRQ(MatMultTransposeAdd(Na->A,Na->w,v2,v3));
178   }
179   PetscFunctionReturn(0);
180 }
181 
182 PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y)
183 {
184   Mat_Normal     *Na = (Mat_Normal*)N->data;
185   Vec            in;
186 
187   PetscFunctionBegin;
188   in = x;
189   if (Na->left) {
190     if (!Na->leftwork) {
191       CHKERRQ(VecDuplicate(Na->left,&Na->leftwork));
192     }
193     CHKERRQ(VecPointwiseMult(Na->leftwork,Na->left,in));
194     in   = Na->leftwork;
195   }
196   CHKERRQ(MatMult(Na->A,in,Na->w));
197   CHKERRQ(MatMultTranspose(Na->A,Na->w,y));
198   if (Na->right) {
199     CHKERRQ(VecPointwiseMult(y,Na->right,y));
200   }
201   CHKERRQ(VecScale(y,Na->scale));
202   PetscFunctionReturn(0);
203 }
204 
205 PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
206 {
207   Mat_Normal     *Na = (Mat_Normal*)N->data;
208   Vec            in;
209 
210   PetscFunctionBegin;
211   in = v1;
212   if (Na->left) {
213     if (!Na->leftwork) {
214       CHKERRQ(VecDuplicate(Na->left,&Na->leftwork));
215     }
216     CHKERRQ(VecPointwiseMult(Na->leftwork,Na->left,in));
217     in   = Na->leftwork;
218   }
219   CHKERRQ(MatMult(Na->A,in,Na->w));
220   CHKERRQ(VecScale(Na->w,Na->scale));
221   if (Na->right) {
222     CHKERRQ(MatMultTranspose(Na->A,Na->w,v3));
223     CHKERRQ(VecPointwiseMult(v3,Na->right,v3));
224     CHKERRQ(VecAXPY(v3,1.0,v2));
225   } else {
226     CHKERRQ(MatMultTransposeAdd(Na->A,Na->w,v2,v3));
227   }
228   PetscFunctionReturn(0);
229 }
230 
231 PetscErrorCode MatDestroy_Normal(Mat N)
232 {
233   Mat_Normal     *Na = (Mat_Normal*)N->data;
234 
235   PetscFunctionBegin;
236   CHKERRQ(MatDestroy(&Na->A));
237   CHKERRQ(VecDestroy(&Na->w));
238   CHKERRQ(VecDestroy(&Na->left));
239   CHKERRQ(VecDestroy(&Na->right));
240   CHKERRQ(VecDestroy(&Na->leftwork));
241   CHKERRQ(VecDestroy(&Na->rightwork));
242   CHKERRQ(PetscFree(N->data));
243   CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatNormalGetMat_C",NULL));
244   CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_seqaij_C",NULL));
245   CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_mpiaij_C",NULL));
246   CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_seqdense_C",NULL));
247   CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_mpidense_C",NULL));
248   CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_dense_C",NULL));
249   PetscFunctionReturn(0);
250 }
251 
252 /*
253       Slow, nonscalable version
254 */
255 PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v)
256 {
257   Mat_Normal        *Na = (Mat_Normal*)N->data;
258   Mat               A   = Na->A;
259   PetscInt          i,j,rstart,rend,nnz;
260   const PetscInt    *cols;
261   PetscScalar       *diag,*work,*values;
262   const PetscScalar *mvalues;
263 
264   PetscFunctionBegin;
265   CHKERRQ(PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work));
266   CHKERRQ(PetscArrayzero(work,A->cmap->N));
267   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
268   for (i=rstart; i<rend; i++) {
269     CHKERRQ(MatGetRow(A,i,&nnz,&cols,&mvalues));
270     for (j=0; j<nnz; j++) {
271       work[cols[j]] += mvalues[j]*mvalues[j];
272     }
273     CHKERRQ(MatRestoreRow(A,i,&nnz,&cols,&mvalues));
274   }
275   CHKERRMPI(MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N)));
276   rstart = N->cmap->rstart;
277   rend   = N->cmap->rend;
278   CHKERRQ(VecGetArray(v,&values));
279   CHKERRQ(PetscArraycpy(values,diag+rstart,rend-rstart));
280   CHKERRQ(VecRestoreArray(v,&values));
281   CHKERRQ(PetscFree2(diag,work));
282   CHKERRQ(VecScale(v,Na->scale));
283   PetscFunctionReturn(0);
284 }
285 
286 PetscErrorCode MatNormalGetMat_Normal(Mat A,Mat *M)
287 {
288   Mat_Normal *Aa = (Mat_Normal*)A->data;
289 
290   PetscFunctionBegin;
291   *M = Aa->A;
292   PetscFunctionReturn(0);
293 }
294 
295 /*@
296       MatNormalGetMat - Gets the Mat object stored inside a MATNORMAL
297 
298    Logically collective on Mat
299 
300    Input Parameter:
301 .   A  - the MATNORMAL matrix
302 
303    Output Parameter:
304 .   M - the matrix object stored inside A
305 
306    Level: intermediate
307 
308 .seealso: MatCreateNormal()
309 
310 @*/
311 PetscErrorCode MatNormalGetMat(Mat A,Mat *M)
312 {
313   PetscFunctionBegin;
314   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
315   PetscValidType(A,1);
316   PetscValidPointer(M,2);
317   CHKERRQ(PetscUseMethod(A,"MatNormalGetMat_C",(Mat,Mat*),(A,M)));
318   PetscFunctionReturn(0);
319 }
320 
321 PetscErrorCode MatConvert_Normal_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
322 {
323   Mat_Normal     *Aa = (Mat_Normal*)A->data;
324   Mat            B;
325   PetscInt       m,n,M,N;
326 
327   PetscFunctionBegin;
328   CHKERRQ(MatGetSize(A,&M,&N));
329   CHKERRQ(MatGetLocalSize(A,&m,&n));
330   if (reuse == MAT_REUSE_MATRIX) {
331     B = *newmat;
332     CHKERRQ(MatProductReplaceMats(Aa->A,Aa->A,NULL,B));
333   } else {
334     CHKERRQ(MatProductCreate(Aa->A,Aa->A,NULL,&B));
335     CHKERRQ(MatProductSetType(B,MATPRODUCT_AtB));
336     CHKERRQ(MatProductSetFromOptions(B));
337     CHKERRQ(MatProductSymbolic(B));
338     CHKERRQ(MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE));
339   }
340   CHKERRQ(MatProductNumeric(B));
341   if (reuse == MAT_INPLACE_MATRIX) {
342     CHKERRQ(MatHeaderReplace(A,&B));
343   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
344   CHKERRQ(MatConvert(*newmat,MATAIJ,MAT_INPLACE_MATRIX,newmat));
345   PetscFunctionReturn(0);
346 }
347 
348 typedef struct {
349   Mat work[2];
350 } Normal_Dense;
351 
352 PetscErrorCode MatProductNumeric_Normal_Dense(Mat C)
353 {
354   Mat            A,B;
355   Normal_Dense   *contents;
356   Mat_Normal     *a;
357   PetscScalar    *array;
358 
359   PetscFunctionBegin;
360   MatCheckProduct(C,3);
361   A = C->product->A;
362   a = (Mat_Normal*)A->data;
363   B = C->product->B;
364   contents = (Normal_Dense*)C->product->data;
365   PetscCheck(contents,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
366   if (a->right) {
367     CHKERRQ(MatCopy(B,C,SAME_NONZERO_PATTERN));
368     CHKERRQ(MatDiagonalScale(C,a->right,NULL));
369   }
370   CHKERRQ(MatProductNumeric(contents->work[0]));
371   CHKERRQ(MatDenseGetArrayWrite(C,&array));
372   CHKERRQ(MatDensePlaceArray(contents->work[1],array));
373   CHKERRQ(MatProductNumeric(contents->work[1]));
374   CHKERRQ(MatDenseRestoreArrayWrite(C,&array));
375   CHKERRQ(MatDenseResetArray(contents->work[1]));
376   CHKERRQ(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE));
377   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
378   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
379   CHKERRQ(MatScale(C,a->scale));
380   PetscFunctionReturn(0);
381 }
382 
383 PetscErrorCode MatNormal_DenseDestroy(void *ctx)
384 {
385   Normal_Dense   *contents = (Normal_Dense*)ctx;
386 
387   PetscFunctionBegin;
388   CHKERRQ(MatDestroy(contents->work));
389   CHKERRQ(MatDestroy(contents->work+1));
390   CHKERRQ(PetscFree(contents));
391   PetscFunctionReturn(0);
392 }
393 
394 PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C)
395 {
396   Mat            A,B;
397   Normal_Dense   *contents = NULL;
398   Mat_Normal     *a;
399   PetscScalar    *array;
400   PetscInt       n,N,m,M;
401 
402   PetscFunctionBegin;
403   MatCheckProduct(C,4);
404   PetscCheck(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
405   A = C->product->A;
406   a = (Mat_Normal*)A->data;
407   PetscCheck(!a->left,PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Not implemented");
408   B = C->product->B;
409   CHKERRQ(MatGetLocalSize(C,&m,&n));
410   CHKERRQ(MatGetSize(C,&M,&N));
411   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
412     CHKERRQ(MatGetLocalSize(B,NULL,&n));
413     CHKERRQ(MatGetSize(B,NULL,&N));
414     CHKERRQ(MatGetLocalSize(A,&m,NULL));
415     CHKERRQ(MatGetSize(A,&M,NULL));
416     CHKERRQ(MatSetSizes(C,m,n,M,N));
417   }
418   CHKERRQ(MatSetType(C,((PetscObject)B)->type_name));
419   CHKERRQ(MatSetUp(C));
420   CHKERRQ(PetscNew(&contents));
421   C->product->data = contents;
422   C->product->destroy = MatNormal_DenseDestroy;
423   if (a->right) {
424     CHKERRQ(MatProductCreate(a->A,C,NULL,contents->work));
425   } else {
426     CHKERRQ(MatProductCreate(a->A,B,NULL,contents->work));
427   }
428   CHKERRQ(MatProductSetType(contents->work[0],MATPRODUCT_AB));
429   CHKERRQ(MatProductSetFromOptions(contents->work[0]));
430   CHKERRQ(MatProductSymbolic(contents->work[0]));
431   CHKERRQ(MatProductCreate(a->A,contents->work[0],NULL,contents->work+1));
432   CHKERRQ(MatProductSetType(contents->work[1],MATPRODUCT_AtB));
433   CHKERRQ(MatProductSetFromOptions(contents->work[1]));
434   CHKERRQ(MatProductSymbolic(contents->work[1]));
435   CHKERRQ(MatDenseGetArrayWrite(C,&array));
436   CHKERRQ(MatSeqDenseSetPreallocation(contents->work[1],array));
437   CHKERRQ(MatMPIDenseSetPreallocation(contents->work[1],array));
438   CHKERRQ(MatDenseRestoreArrayWrite(C,&array));
439   C->ops->productnumeric = MatProductNumeric_Normal_Dense;
440   PetscFunctionReturn(0);
441 }
442 
443 PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C)
444 {
445   PetscFunctionBegin;
446   C->ops->productsymbolic = MatProductSymbolic_Normal_Dense;
447   PetscFunctionReturn(0);
448 }
449 
450 PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C)
451 {
452   Mat_Product    *product = C->product;
453 
454   PetscFunctionBegin;
455   if (product->type == MATPRODUCT_AB) {
456     CHKERRQ(MatProductSetFromOptions_Normal_Dense_AB(C));
457   }
458   PetscFunctionReturn(0);
459 }
460 
461 /*@
462       MatCreateNormal - Creates a new matrix object that behaves like A'*A.
463 
464    Collective on Mat
465 
466    Input Parameter:
467 .   A  - the (possibly rectangular) matrix
468 
469    Output Parameter:
470 .   N - the matrix that represents A'*A
471 
472    Level: intermediate
473 
474    Notes:
475     The product A'*A is NOT actually formed! Rather the new matrix
476           object performs the matrix-vector product by first multiplying by
477           A and then A'
478 @*/
479 PetscErrorCode  MatCreateNormal(Mat A,Mat *N)
480 {
481   PetscInt       n,nn;
482   Mat_Normal     *Na;
483   VecType        vtype;
484 
485   PetscFunctionBegin;
486   CHKERRQ(MatGetSize(A,NULL,&nn));
487   CHKERRQ(MatGetLocalSize(A,NULL,&n));
488   CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),N));
489   CHKERRQ(MatSetSizes(*N,n,n,nn,nn));
490   CHKERRQ(PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL));
491   CHKERRQ(PetscLayoutReference(A->cmap,&(*N)->rmap));
492   CHKERRQ(PetscLayoutReference(A->cmap,&(*N)->cmap));
493 
494   CHKERRQ(PetscNewLog(*N,&Na));
495   (*N)->data = (void*) Na;
496   CHKERRQ(PetscObjectReference((PetscObject)A));
497   Na->A      = A;
498   Na->scale  = 1.0;
499 
500   CHKERRQ(MatCreateVecs(A,NULL,&Na->w));
501 
502   (*N)->ops->destroy           = MatDestroy_Normal;
503   (*N)->ops->mult              = MatMult_Normal;
504   (*N)->ops->multtranspose     = MatMultTranspose_Normal;
505   (*N)->ops->multtransposeadd  = MatMultTransposeAdd_Normal;
506   (*N)->ops->multadd           = MatMultAdd_Normal;
507   (*N)->ops->getdiagonal       = MatGetDiagonal_Normal;
508   (*N)->ops->scale             = MatScale_Normal;
509   (*N)->ops->diagonalscale     = MatDiagonalScale_Normal;
510   (*N)->ops->increaseoverlap   = MatIncreaseOverlap_Normal;
511   (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal;
512   (*N)->ops->permute           = MatPermute_Normal;
513   (*N)->ops->duplicate         = MatDuplicate_Normal;
514   (*N)->ops->copy              = MatCopy_Normal;
515   (*N)->assembled              = PETSC_TRUE;
516   (*N)->preallocated           = PETSC_TRUE;
517 
518   CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatNormalGetMat_C",MatNormalGetMat_Normal));
519   CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_seqaij_C",MatConvert_Normal_AIJ));
520   CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_mpiaij_C",MatConvert_Normal_AIJ));
521   CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_seqdense_C",MatProductSetFromOptions_Normal_Dense));
522   CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_mpidense_C",MatProductSetFromOptions_Normal_Dense));
523   CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_dense_C",MatProductSetFromOptions_Normal_Dense));
524   CHKERRQ(MatSetOption(*N,MAT_SYMMETRIC,PETSC_TRUE));
525   CHKERRQ(MatGetVecType(A,&vtype));
526   CHKERRQ(MatSetVecType(*N,vtype));
527 #if defined(PETSC_HAVE_DEVICE)
528   CHKERRQ(MatBindToCPU(*N,A->boundtocpu));
529 #endif
530   PetscFunctionReturn(0);
531 }
532