1 2 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3 4 /*@ 5 MatAXPY - Computes Y = a*X + Y. 6 7 Logically Collective on Mat 8 9 Input Parameters: 10 + a - the scalar multiplier 11 . X - the first matrix 12 . Y - the second matrix 13 - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN 14 or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's) 15 16 Level: intermediate 17 18 .keywords: matrix, add 19 20 .seealso: MatAYPX() 21 @*/ 22 PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str) 23 { 24 PetscErrorCode ierr; 25 PetscInt m1,m2,n1,n2; 26 PetscBool sametype; 27 28 PetscFunctionBegin; 29 PetscValidHeaderSpecific(X,MAT_CLASSID,3); 30 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 31 PetscValidLogicalCollectiveScalar(Y,a,2); 32 ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr); 33 ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr); 34 if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrix add: %D %D %D %D",m1,m2,n1,n2); 35 36 ierr = PetscStrcmp(((PetscObject)X)->type_name,((PetscObject)Y)->type_name,&sametype);CHKERRQ(ierr); 37 ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 38 if (Y->ops->axpy && sametype) { 39 ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr); 40 } else { 41 ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 42 } 43 ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 44 #if defined(PETSC_HAVE_CUSP) 45 if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) { 46 Y->valid_GPU_matrix = PETSC_CUSP_CPU; 47 } 48 #elif defined(PETSC_HAVE_VIENNACL) 49 if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) { 50 Y->valid_GPU_matrix = PETSC_VIENNACL_CPU; 51 } 52 #elif defined(PETSC_HAVE_VECCUDA) 53 if (Y->valid_GPU_matrix != PETSC_CUDA_UNALLOCATED) { 54 Y->valid_GPU_matrix = PETSC_CUDA_CPU; 55 } 56 #endif 57 PetscFunctionReturn(0); 58 } 59 60 PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 61 { 62 PetscInt i,start,end,j,ncols,m,n; 63 PetscErrorCode ierr; 64 const PetscInt *row; 65 PetscScalar *val; 66 const PetscScalar *vals; 67 68 PetscFunctionBegin; 69 ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 70 ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 71 if (a == 1.0) { 72 for (i = start; i < end; i++) { 73 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 74 ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 75 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 76 } 77 } else { 78 ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr); 79 for (i=start; i<end; i++) { 80 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 81 for (j=0; j<ncols; j++) { 82 val[j] = a*vals[j]; 83 } 84 ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 85 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 86 } 87 ierr = PetscFree(val);CHKERRQ(ierr); 88 } 89 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 90 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 91 PetscFunctionReturn(0); 92 } 93 94 PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 95 { 96 PetscInt i,start,end,j,ncols,m,n; 97 PetscErrorCode ierr; 98 const PetscInt *row; 99 PetscScalar *val; 100 const PetscScalar *vals; 101 102 PetscFunctionBegin; 103 ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 104 ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 105 if (a == 1.0) { 106 for (i = start; i < end; i++) { 107 ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 108 ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 109 ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 110 111 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 112 ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 113 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 114 } 115 } else { 116 ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr); 117 for (i=start; i<end; i++) { 118 ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 119 ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 120 ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 121 122 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 123 for (j=0; j<ncols; j++) { 124 val[j] = a*vals[j]; 125 } 126 ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 127 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 128 } 129 ierr = PetscFree(val);CHKERRQ(ierr); 130 } 131 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 132 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133 PetscFunctionReturn(0); 134 } 135 136 /*@ 137 MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 138 139 Neighbor-wise Collective on Mat 140 141 Input Parameters: 142 + Y - the matrices 143 - a - the PetscScalar 144 145 Level: intermediate 146 147 Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 148 fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an 149 entry). 150 151 To form Y = Y + diag(V) use MatDiagonalSet() 152 153 Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is 154 preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix. 155 156 .keywords: matrix, add, shift 157 158 .seealso: MatDiagonalSet(), MatScale(), MatDiagonalScale() 159 @*/ 160 PetscErrorCode MatShift(Mat Y,PetscScalar a) 161 { 162 PetscErrorCode ierr; 163 164 PetscFunctionBegin; 165 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 166 if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 167 if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 168 MatCheckPreallocated(Y,1); 169 170 if (Y->ops->shift) { 171 ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 172 } else { 173 ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 174 } 175 176 #if defined(PETSC_HAVE_CUSP) 177 if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) { 178 Y->valid_GPU_matrix = PETSC_CUSP_CPU; 179 } 180 #elif defined(PETSC_HAVE_VIENNACL) 181 if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) { 182 Y->valid_GPU_matrix = PETSC_VIENNACL_CPU; 183 } 184 #elif defined(PETSC_HAVE_VECCUDA) 185 if (Y->valid_GPU_matrix != PETSC_CUDA_UNALLOCATED) { 186 Y->valid_GPU_matrix = PETSC_CUDA_CPU; 187 } 188 #endif 189 PetscFunctionReturn(0); 190 } 191 192 PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 193 { 194 PetscErrorCode ierr; 195 PetscInt i,start,end; 196 PetscScalar *v; 197 198 PetscFunctionBegin; 199 ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 200 ierr = VecGetArray(D,&v);CHKERRQ(ierr); 201 for (i=start; i<end; i++) { 202 ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 203 } 204 ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 205 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 206 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 207 PetscFunctionReturn(0); 208 } 209 210 /*@ 211 MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 212 that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 213 INSERT_VALUES. 214 215 Input Parameters: 216 + Y - the input matrix 217 . D - the diagonal matrix, represented as a vector 218 - i - INSERT_VALUES or ADD_VALUES 219 220 Neighbor-wise Collective on Mat and Vec 221 222 Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 223 fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an 224 entry). 225 226 Level: intermediate 227 228 .keywords: matrix, add, shift, diagonal 229 230 .seealso: MatShift(), MatScale(), MatDiagonalScale() 231 @*/ 232 PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 233 { 234 PetscErrorCode ierr; 235 PetscInt matlocal,veclocal; 236 237 PetscFunctionBegin; 238 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 239 PetscValidHeaderSpecific(D,VEC_CLASSID,2); 240 ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 241 ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 242 if (matlocal != veclocal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number local rows of matrix %D does not match that of vector for diagonal %D",matlocal,veclocal); 243 if (Y->ops->diagonalset) { 244 ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 245 } else { 246 ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 247 } 248 PetscFunctionReturn(0); 249 } 250 251 /*@ 252 MatAYPX - Computes Y = a*Y + X. 253 254 Logically on Mat 255 256 Input Parameters: 257 + a - the PetscScalar multiplier 258 . Y - the first matrix 259 . X - the second matrix 260 - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 261 262 Level: intermediate 263 264 .keywords: matrix, add 265 266 .seealso: MatAXPY() 267 @*/ 268 PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 269 { 270 PetscScalar one = 1.0; 271 PetscErrorCode ierr; 272 PetscInt mX,mY,nX,nY; 273 274 PetscFunctionBegin; 275 PetscValidHeaderSpecific(X,MAT_CLASSID,3); 276 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 277 PetscValidLogicalCollectiveScalar(Y,a,2); 278 ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 279 ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 280 if (mX != mY || nX != nY) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrices: %D %D first %D %D second",mX,mY,nX,nY); 281 282 ierr = MatScale(Y,a);CHKERRQ(ierr); 283 ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 /*@ 288 MatComputeExplicitOperator - Computes the explicit matrix 289 290 Collective on Mat 291 292 Input Parameter: 293 . inmat - the matrix 294 295 Output Parameter: 296 . mat - the explict operator 297 298 Notes: 299 This computation is done by applying the operators to columns of the 300 identity matrix. 301 302 Currently, this routine uses a dense matrix format when 1 processor 303 is used and a sparse format otherwise. This routine is costly in general, 304 and is recommended for use only with relatively small systems. 305 306 Level: advanced 307 308 .keywords: Mat, compute, explicit, operator 309 @*/ 310 PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 311 { 312 Vec in,out; 313 PetscErrorCode ierr; 314 PetscInt i,m,n,M,N,*rows,start,end; 315 MPI_Comm comm; 316 PetscScalar *array,zero = 0.0,one = 1.0; 317 PetscMPIInt size; 318 319 PetscFunctionBegin; 320 PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 321 PetscValidPointer(mat,2); 322 323 ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 324 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 325 326 ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 327 ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 328 ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 329 ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 330 ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 331 ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 332 for (i=0; i<m; i++) rows[i] = start + i; 333 334 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 335 ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 336 if (size == 1) { 337 ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 338 ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 339 } else { 340 ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 341 ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 342 } 343 344 for (i=0; i<N; i++) { 345 346 ierr = VecSet(in,zero);CHKERRQ(ierr); 347 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 348 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 349 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 350 351 ierr = MatMult(inmat,in,out);CHKERRQ(ierr); 352 353 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 354 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 355 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 356 357 } 358 ierr = PetscFree(rows);CHKERRQ(ierr); 359 ierr = VecDestroy(&out);CHKERRQ(ierr); 360 ierr = VecDestroy(&in);CHKERRQ(ierr); 361 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 362 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 363 PetscFunctionReturn(0); 364 } 365 366 /*@ 367 MatComputeExplicitOperatorTranspose - Computes the explicit matrix representation of 368 a give matrix that can apply MatMultTranspose() 369 370 Collective on Mat 371 372 Input Parameter: 373 . inmat - the matrix 374 375 Output Parameter: 376 . mat - the explict operator transposed 377 378 Notes: 379 This computation is done by applying the transpose of the operator to columns of the 380 identity matrix. 381 382 Currently, this routine uses a dense matrix format when 1 processor 383 is used and a sparse format otherwise. This routine is costly in general, 384 and is recommended for use only with relatively small systems. 385 386 Level: advanced 387 388 .keywords: Mat, compute, explicit, operator 389 @*/ 390 PetscErrorCode MatComputeExplicitOperatorTranspose(Mat inmat,Mat *mat) 391 { 392 Vec in,out; 393 PetscErrorCode ierr; 394 PetscInt i,m,n,M,N,*rows,start,end; 395 MPI_Comm comm; 396 PetscScalar *array,zero = 0.0,one = 1.0; 397 PetscMPIInt size; 398 399 PetscFunctionBegin; 400 PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 401 PetscValidPointer(mat,2); 402 403 ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 404 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 405 406 ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 407 ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 408 ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 409 ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 410 ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 411 ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 412 for (i=0; i<m; i++) rows[i] = start + i; 413 414 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 415 ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 416 if (size == 1) { 417 ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 418 ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 419 } else { 420 ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 421 ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 422 } 423 424 for (i=0; i<N; i++) { 425 426 ierr = VecSet(in,zero);CHKERRQ(ierr); 427 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 428 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 429 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 430 431 ierr = MatMultTranspose(inmat,in,out);CHKERRQ(ierr); 432 433 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 434 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 435 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 436 437 } 438 ierr = PetscFree(rows);CHKERRQ(ierr); 439 ierr = VecDestroy(&out);CHKERRQ(ierr); 440 ierr = VecDestroy(&in);CHKERRQ(ierr); 441 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 442 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 443 PetscFunctionReturn(0); 444 } 445 446 /*@ 447 MatChop - Set all values in the matrix less than the tolerance to zero 448 449 Input Parameters: 450 + A - The matrix 451 - tol - The zero tolerance 452 453 Output Parameters: 454 . A - The chopped matrix 455 456 Level: intermediate 457 458 .seealso: MatCreate(), MatZeroEntries() 459 @*/ 460 PetscErrorCode MatChop(Mat A, PetscReal tol) 461 { 462 PetscScalar *newVals; 463 PetscInt *newCols; 464 PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 465 PetscErrorCode ierr; 466 467 PetscFunctionBegin; 468 ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 469 for (r = rStart; r < rEnd; ++r) { 470 PetscInt ncols; 471 472 ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 473 colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 474 ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 475 } 476 numRows = rEnd - rStart; 477 ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 478 ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 479 for (r = rStart; r < rStart+maxRows; ++r) { 480 const PetscScalar *vals; 481 const PetscInt *cols; 482 PetscInt ncols, newcols, c; 483 484 if (r < rEnd) { 485 ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 486 for (c = 0; c < ncols; ++c) { 487 newCols[c] = cols[c]; 488 newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 489 } 490 newcols = ncols; 491 ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 492 ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 493 } 494 ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 495 ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 496 } 497 ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 498 PetscFunctionReturn(0); 499 } 500