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 Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is 152 preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix. 153 154 .keywords: matrix, add, shift 155 156 .seealso: MatDiagonalSet() 157 @*/ 158 PetscErrorCode MatShift(Mat Y,PetscScalar a) 159 { 160 PetscErrorCode ierr; 161 162 PetscFunctionBegin; 163 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 164 if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 165 if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 166 MatCheckPreallocated(Y,1); 167 168 if (Y->ops->shift) { 169 ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 170 } else { 171 ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 172 } 173 174 #if defined(PETSC_HAVE_CUSP) 175 if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) { 176 Y->valid_GPU_matrix = PETSC_CUSP_CPU; 177 } 178 #elif defined(PETSC_HAVE_VIENNACL) 179 if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) { 180 Y->valid_GPU_matrix = PETSC_VIENNACL_CPU; 181 } 182 #elif defined(PETSC_HAVE_VECCUDA) 183 if (Y->valid_GPU_matrix != PETSC_CUDA_UNALLOCATED) { 184 Y->valid_GPU_matrix = PETSC_CUDA_CPU; 185 } 186 #endif 187 PetscFunctionReturn(0); 188 } 189 190 PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 191 { 192 PetscErrorCode ierr; 193 PetscInt i,start,end; 194 PetscScalar *v; 195 196 PetscFunctionBegin; 197 ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 198 ierr = VecGetArray(D,&v);CHKERRQ(ierr); 199 for (i=start; i<end; i++) { 200 ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 201 } 202 ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 203 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 204 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 205 PetscFunctionReturn(0); 206 } 207 208 /*@ 209 MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 210 that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 211 INSERT_VALUES. 212 213 Input Parameters: 214 + Y - the input matrix 215 . D - the diagonal matrix, represented as a vector 216 - i - INSERT_VALUES or ADD_VALUES 217 218 Neighbor-wise Collective on Mat and Vec 219 220 Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 221 fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an 222 entry). 223 224 Level: intermediate 225 226 .keywords: matrix, add, shift, diagonal 227 228 .seealso: MatShift() 229 @*/ 230 PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 231 { 232 PetscErrorCode ierr; 233 PetscInt matlocal,veclocal; 234 235 PetscFunctionBegin; 236 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 237 PetscValidHeaderSpecific(D,VEC_CLASSID,2); 238 ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 239 ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 240 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); 241 if (Y->ops->diagonalset) { 242 ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 243 } else { 244 ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 245 } 246 PetscFunctionReturn(0); 247 } 248 249 /*@ 250 MatAYPX - Computes Y = a*Y + X. 251 252 Logically on Mat 253 254 Input Parameters: 255 + a - the PetscScalar multiplier 256 . Y - the first matrix 257 . X - the second matrix 258 - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 259 260 Level: intermediate 261 262 .keywords: matrix, add 263 264 .seealso: MatAXPY() 265 @*/ 266 PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 267 { 268 PetscScalar one = 1.0; 269 PetscErrorCode ierr; 270 PetscInt mX,mY,nX,nY; 271 272 PetscFunctionBegin; 273 PetscValidHeaderSpecific(X,MAT_CLASSID,3); 274 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 275 PetscValidLogicalCollectiveScalar(Y,a,2); 276 ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 277 ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 278 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); 279 280 ierr = MatScale(Y,a);CHKERRQ(ierr); 281 ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 282 PetscFunctionReturn(0); 283 } 284 285 /*@ 286 MatComputeExplicitOperator - Computes the explicit matrix 287 288 Collective on Mat 289 290 Input Parameter: 291 . inmat - the matrix 292 293 Output Parameter: 294 . mat - the explict operator 295 296 Notes: 297 This computation is done by applying the operators to columns of the 298 identity matrix. 299 300 Currently, this routine uses a dense matrix format when 1 processor 301 is used and a sparse format otherwise. This routine is costly in general, 302 and is recommended for use only with relatively small systems. 303 304 Level: advanced 305 306 .keywords: Mat, compute, explicit, operator 307 @*/ 308 PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 309 { 310 Vec in,out; 311 PetscErrorCode ierr; 312 PetscInt i,m,n,M,N,*rows,start,end; 313 MPI_Comm comm; 314 PetscScalar *array,zero = 0.0,one = 1.0; 315 PetscMPIInt size; 316 317 PetscFunctionBegin; 318 PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 319 PetscValidPointer(mat,2); 320 321 ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 322 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 323 324 ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 325 ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 326 ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 327 ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 328 ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 329 ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 330 for (i=0; i<m; i++) rows[i] = start + i; 331 332 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 333 ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 334 if (size == 1) { 335 ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 336 ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 337 } else { 338 ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 339 ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 340 } 341 342 for (i=0; i<N; i++) { 343 344 ierr = VecSet(in,zero);CHKERRQ(ierr); 345 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 346 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 347 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 348 349 ierr = MatMult(inmat,in,out);CHKERRQ(ierr); 350 351 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 352 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 353 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 354 355 } 356 ierr = PetscFree(rows);CHKERRQ(ierr); 357 ierr = VecDestroy(&out);CHKERRQ(ierr); 358 ierr = VecDestroy(&in);CHKERRQ(ierr); 359 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 360 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 361 PetscFunctionReturn(0); 362 } 363 364 /*@ 365 MatComputeExplicitOperatorTranspose - Computes the explicit matrix representation of 366 a give matrix that can apply MatMultTranspose() 367 368 Collective on Mat 369 370 Input Parameter: 371 . inmat - the matrix 372 373 Output Parameter: 374 . mat - the explict operator transposed 375 376 Notes: 377 This computation is done by applying the transpose of the operator to columns of the 378 identity matrix. 379 380 Currently, this routine uses a dense matrix format when 1 processor 381 is used and a sparse format otherwise. This routine is costly in general, 382 and is recommended for use only with relatively small systems. 383 384 Level: advanced 385 386 .keywords: Mat, compute, explicit, operator 387 @*/ 388 PetscErrorCode MatComputeExplicitOperatorTranspose(Mat inmat,Mat *mat) 389 { 390 Vec in,out; 391 PetscErrorCode ierr; 392 PetscInt i,m,n,M,N,*rows,start,end; 393 MPI_Comm comm; 394 PetscScalar *array,zero = 0.0,one = 1.0; 395 PetscMPIInt size; 396 397 PetscFunctionBegin; 398 PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 399 PetscValidPointer(mat,2); 400 401 ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 402 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 403 404 ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 405 ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 406 ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 407 ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 408 ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 409 ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 410 for (i=0; i<m; i++) rows[i] = start + i; 411 412 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 413 ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 414 if (size == 1) { 415 ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 416 ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 417 } else { 418 ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 419 ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 420 } 421 422 for (i=0; i<N; i++) { 423 424 ierr = VecSet(in,zero);CHKERRQ(ierr); 425 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 426 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 427 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 428 429 ierr = MatMultTranspose(inmat,in,out);CHKERRQ(ierr); 430 431 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 432 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 433 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 434 435 } 436 ierr = PetscFree(rows);CHKERRQ(ierr); 437 ierr = VecDestroy(&out);CHKERRQ(ierr); 438 ierr = VecDestroy(&in);CHKERRQ(ierr); 439 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 440 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 441 PetscFunctionReturn(0); 442 } 443 444 /*@ 445 MatChop - Set all values in the matrix less than the tolerance to zero 446 447 Input Parameters: 448 + A - The matrix 449 - tol - The zero tolerance 450 451 Output Parameters: 452 . A - The chopped matrix 453 454 Level: intermediate 455 456 .seealso: MatCreate(), MatZeroEntries() 457 @*/ 458 PetscErrorCode MatChop(Mat A, PetscReal tol) 459 { 460 PetscScalar *newVals; 461 PetscInt *newCols; 462 PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 463 PetscErrorCode ierr; 464 465 PetscFunctionBegin; 466 ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 467 for (r = rStart; r < rEnd; ++r) { 468 PetscInt ncols; 469 470 ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 471 colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 472 ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 473 } 474 numRows = rEnd - rStart; 475 ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 476 ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 477 for (r = rStart; r < rStart+maxRows; ++r) { 478 const PetscScalar *vals; 479 const PetscInt *cols; 480 PetscInt ncols, newcols, c; 481 482 if (r < rEnd) { 483 ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 484 for (c = 0; c < ncols; ++c) { 485 newCols[c] = cols[c]; 486 newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 487 } 488 newcols = ncols; 489 ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 490 ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 491 } 492 ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 493 ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 494 } 495 ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 496 PetscFunctionReturn(0); 497 } 498