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 PetscScalar shift, scale; 32 PetscInt M; 33 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 } 56 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; 63 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 } 79 80 static PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B) 81 { 82 Mat_Normal *a; 83 Mat C; 84 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 } 93 94 static PetscErrorCode MatCopy_Normal(Mat A, Mat B, MatStructure str) 95 { 96 Mat_Normal *a, *b; 97 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 } 104 105 static PetscErrorCode MatMult_Normal(Mat N, Vec x, Vec y) 106 { 107 Mat_Normal *Na; 108 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 } 115 116 static PetscErrorCode MatDestroy_Normal(Mat N) 117 { 118 Mat_Normal *Na; 119 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 } 137 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 *work, *values; 148 const PetscScalar *mvalues; 149 150 PetscFunctionBegin; 151 PetscCall(MatShellGetContext(N, &Na)); 152 A = Na->A; 153 PetscCall(PetscMalloc1(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 PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, work, 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, work + rstart, rend - rstart)); 166 PetscCall(VecRestoreArray(v, &values)); 167 PetscCall(PetscFree(work)); 168 PetscFunctionReturn(PETSC_SUCCESS); 169 } 170 171 static PetscErrorCode MatGetDiagonalBlock_Normal(Mat N, Mat *D) 172 { 173 Mat_Normal *Na; 174 Mat M, A; 175 176 PetscFunctionBegin; 177 PetscCall(MatShellGetContext(N, &Na)); 178 A = Na->A; 179 PetscCall(MatGetDiagonalBlock(A, &M)); 180 PetscCall(MatCreateNormal(M, &Na->D)); 181 *D = Na->D; 182 PetscFunctionReturn(PETSC_SUCCESS); 183 } 184 185 static PetscErrorCode MatNormalGetMat_Normal(Mat A, Mat *M) 186 { 187 Mat_Normal *Aa; 188 189 PetscFunctionBegin; 190 PetscCall(MatShellGetContext(A, &Aa)); 191 *M = Aa->A; 192 PetscFunctionReturn(PETSC_SUCCESS); 193 } 194 195 /*@ 196 MatNormalGetMat - Gets the `Mat` object stored inside a `MATNORMAL` 197 198 Logically Collective 199 200 Input Parameter: 201 . A - the `MATNORMAL` matrix 202 203 Output Parameter: 204 . M - the matrix object stored inside `A` 205 206 Level: intermediate 207 208 .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MATNORMALHERMITIAN`, `MatCreateNormal()` 209 @*/ 210 PetscErrorCode MatNormalGetMat(Mat A, Mat *M) 211 { 212 PetscFunctionBegin; 213 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 214 PetscValidType(A, 1); 215 PetscAssertPointer(M, 2); 216 PetscUseMethod(A, "MatNormalGetMat_C", (Mat, Mat *), (A, M)); 217 PetscFunctionReturn(PETSC_SUCCESS); 218 } 219 220 static PetscErrorCode MatConvert_Normal_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) 221 { 222 Mat_Normal *Aa; 223 Mat B; 224 PetscInt m, n, M, N; 225 226 PetscFunctionBegin; 227 PetscCall(MatShellGetContext(A, &Aa)); 228 PetscCall(MatGetSize(A, &M, &N)); 229 PetscCall(MatGetLocalSize(A, &m, &n)); 230 if (reuse == MAT_REUSE_MATRIX) { 231 B = *newmat; 232 PetscCall(MatProductReplaceMats(Aa->A, Aa->A, NULL, B)); 233 } else { 234 PetscCall(MatProductCreate(Aa->A, Aa->A, NULL, &B)); 235 PetscCall(MatProductSetType(B, MATPRODUCT_AtB)); 236 PetscCall(MatProductSetFromOptions(B)); 237 PetscCall(MatProductSymbolic(B)); 238 PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE)); 239 } 240 PetscCall(MatProductNumeric(B)); 241 if (reuse == MAT_INPLACE_MATRIX) { 242 PetscCall(MatHeaderReplace(A, &B)); 243 } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 244 PetscCall(MatConvert(*newmat, MATAIJ, MAT_INPLACE_MATRIX, newmat)); 245 PetscFunctionReturn(PETSC_SUCCESS); 246 } 247 248 #if defined(PETSC_HAVE_HYPRE) 249 static PetscErrorCode MatConvert_Normal_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) 250 { 251 PetscFunctionBegin; 252 if (reuse == MAT_INITIAL_MATRIX) { 253 PetscCall(MatConvert(A, MATAIJ, reuse, B)); 254 PetscCall(MatConvert(*B, type, MAT_INPLACE_MATRIX, B)); 255 } else PetscCall(MatConvert_Basic(A, type, reuse, B)); /* fall back to basic convert */ 256 PetscFunctionReturn(PETSC_SUCCESS); 257 } 258 #endif 259 260 typedef struct { 261 Mat work[2]; 262 } Normal_Dense; 263 264 static PetscErrorCode MatProductNumeric_Normal_Dense(Mat C) 265 { 266 Mat A, B; 267 Normal_Dense *contents; 268 Mat_Normal *a; 269 Vec right; 270 PetscScalar *array, scale; 271 272 PetscFunctionBegin; 273 MatCheckProduct(C, 1); 274 A = C->product->A; 275 B = C->product->B; 276 PetscCall(MatShellGetContext(A, &a)); 277 contents = (Normal_Dense *)C->product->data; 278 PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty"); 279 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)); 280 if (right) { 281 PetscCall(MatCopy(B, C, SAME_NONZERO_PATTERN)); 282 PetscCall(MatDiagonalScale(C, right, NULL)); 283 } 284 PetscCall(MatProductNumeric(contents->work[0])); 285 PetscCall(MatDenseGetArrayWrite(C, &array)); 286 PetscCall(MatDensePlaceArray(contents->work[1], array)); 287 PetscCall(MatProductNumeric(contents->work[1])); 288 PetscCall(MatDenseRestoreArrayWrite(C, &array)); 289 PetscCall(MatDenseResetArray(contents->work[1])); 290 PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 291 PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 292 PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 293 PetscCall(MatScale(C, scale)); 294 PetscFunctionReturn(PETSC_SUCCESS); 295 } 296 297 static PetscErrorCode MatNormal_DenseDestroy(void *ctx) 298 { 299 Normal_Dense *contents = (Normal_Dense *)ctx; 300 301 PetscFunctionBegin; 302 PetscCall(MatDestroy(contents->work)); 303 PetscCall(MatDestroy(contents->work + 1)); 304 PetscCall(PetscFree(contents)); 305 PetscFunctionReturn(PETSC_SUCCESS); 306 } 307 308 static PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C) 309 { 310 Mat A, B; 311 Normal_Dense *contents = NULL; 312 Mat_Normal *a; 313 Vec right; 314 PetscScalar *array, scale; 315 PetscInt n, N, m, M; 316 317 PetscFunctionBegin; 318 MatCheckProduct(C, 1); 319 PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty"); 320 A = C->product->A; 321 B = C->product->B; 322 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)); 323 PetscCall(MatShellGetContext(A, &a)); 324 PetscCall(MatGetLocalSize(C, &m, &n)); 325 PetscCall(MatGetSize(C, &M, &N)); 326 if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) { 327 PetscCall(MatGetLocalSize(B, NULL, &n)); 328 PetscCall(MatGetSize(B, NULL, &N)); 329 PetscCall(MatGetLocalSize(A, &m, NULL)); 330 PetscCall(MatGetSize(A, &M, NULL)); 331 PetscCall(MatSetSizes(C, m, n, M, N)); 332 } 333 PetscCall(MatSetType(C, ((PetscObject)B)->type_name)); 334 PetscCall(MatSetUp(C)); 335 PetscCall(PetscNew(&contents)); 336 C->product->data = contents; 337 C->product->destroy = MatNormal_DenseDestroy; 338 if (right) PetscCall(MatProductCreate(a->A, C, NULL, contents->work)); 339 else PetscCall(MatProductCreate(a->A, B, NULL, contents->work)); 340 PetscCall(MatProductSetType(contents->work[0], MATPRODUCT_AB)); 341 PetscCall(MatProductSetFromOptions(contents->work[0])); 342 PetscCall(MatProductSymbolic(contents->work[0])); 343 PetscCall(MatProductCreate(a->A, contents->work[0], NULL, contents->work + 1)); 344 PetscCall(MatProductSetType(contents->work[1], MATPRODUCT_AtB)); 345 PetscCall(MatProductSetFromOptions(contents->work[1])); 346 PetscCall(MatProductSymbolic(contents->work[1])); 347 PetscCall(MatDenseGetArrayWrite(C, &array)); 348 PetscCall(MatSeqDenseSetPreallocation(contents->work[1], array)); 349 PetscCall(MatMPIDenseSetPreallocation(contents->work[1], array)); 350 PetscCall(MatDenseRestoreArrayWrite(C, &array)); 351 C->ops->productnumeric = MatProductNumeric_Normal_Dense; 352 PetscFunctionReturn(PETSC_SUCCESS); 353 } 354 355 static PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C) 356 { 357 Mat_Product *product = C->product; 358 359 PetscFunctionBegin; 360 if (product->type == MATPRODUCT_AB) C->ops->productsymbolic = MatProductSymbolic_Normal_Dense; 361 PetscFunctionReturn(PETSC_SUCCESS); 362 } 363 364 /*MC 365 MATNORMAL - a matrix that behaves like A'*A for `MatMult()` while only containing A 366 367 Level: intermediate 368 369 Developer Notes: 370 This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code 371 372 Users can not call `MatShellSetOperation()` operations on this class, there is some error checking for that incorrect usage 373 374 .seealso: [](ch_matrices), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()` 375 M*/ 376 377 /*@ 378 MatCreateNormal - Creates a new `MATNORMAL` matrix object that behaves like A'*A. 379 380 Collective 381 382 Input Parameter: 383 . A - the (possibly rectangular) matrix 384 385 Output Parameter: 386 . N - the matrix that represents A'*A 387 388 Level: intermediate 389 390 Notes: 391 The product A'*A is NOT actually formed! Rather the new matrix 392 object performs the matrix-vector product, `MatMult()`, by first multiplying by 393 A and then A' 394 395 .seealso: [](ch_matrices), `Mat`, `MATNORMAL`, `MatMult()`, `MatNormalGetMat()`, `MATNORMALHERMITIAN`, `MatCreateNormalHermitian()` 396 @*/ 397 PetscErrorCode MatCreateNormal(Mat A, Mat *N) 398 { 399 Mat_Normal *Na; 400 VecType vtype; 401 402 PetscFunctionBegin; 403 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), N)); 404 PetscCall(PetscLayoutReference(A->cmap, &(*N)->rmap)); 405 PetscCall(PetscLayoutReference(A->cmap, &(*N)->cmap)); 406 PetscCall(MatSetType(*N, MATSHELL)); 407 PetscCall(PetscNew(&Na)); 408 PetscCall(MatShellSetContext(*N, Na)); 409 PetscCall(PetscObjectReference((PetscObject)A)); 410 Na->A = A; 411 PetscCall(MatCreateVecs(A, NULL, &Na->w)); 412 413 PetscCall(MatSetBlockSizes(*N, PetscAbs(A->cmap->bs), PetscAbs(A->rmap->bs))); 414 PetscCall(MatShellSetOperation(*N, MATOP_DESTROY, (void (*)(void))MatDestroy_Normal)); 415 PetscCall(MatShellSetOperation(*N, MATOP_MULT, (void (*)(void))MatMult_Normal)); 416 PetscCall(MatShellSetOperation(*N, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_Normal)); 417 PetscCall(MatShellSetOperation(*N, MATOP_DUPLICATE, (void (*)(void))MatDuplicate_Normal)); 418 PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Normal)); 419 PetscCall(MatShellSetOperation(*N, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Normal)); 420 PetscCall(MatShellSetOperation(*N, MATOP_COPY, (void (*)(void))MatCopy_Normal)); 421 (*N)->ops->increaseoverlap = MatIncreaseOverlap_Normal; 422 (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal; 423 (*N)->ops->permute = MatPermute_Normal; 424 425 PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatNormalGetMat_C", MatNormalGetMat_Normal)); 426 PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_seqaij_C", MatConvert_Normal_AIJ)); 427 PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_mpiaij_C", MatConvert_Normal_AIJ)); 428 #if defined(PETSC_HAVE_HYPRE) 429 PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatConvert_normal_hypre_C", MatConvert_Normal_HYPRE)); 430 #endif 431 PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_seqdense_C", MatProductSetFromOptions_Normal_Dense)); 432 PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatProductSetFromOptions_normal_mpidense_C", MatProductSetFromOptions_Normal_Dense)); 433 PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContext_C", MatShellSetContext_Immutable)); 434 PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable)); 435 PetscCall(PetscObjectComposeFunction((PetscObject)*N, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable)); 436 PetscCall(MatSetOption(*N, MAT_SYMMETRIC, PETSC_TRUE)); 437 PetscCall(MatGetVecType(A, &vtype)); 438 PetscCall(MatSetVecType(*N, vtype)); 439 #if defined(PETSC_HAVE_DEVICE) 440 PetscCall(MatBindToCPU(*N, A->boundtocpu)); 441 #endif 442 PetscCall(MatSetUp(*N)); 443 PetscCall(PetscObjectChangeTypeName((PetscObject)*N, MATNORMAL)); 444 PetscFunctionReturn(PETSC_SUCCESS); 445 } 446