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