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