1 2 /* 3 spbas_cholesky_row_alloc: 4 in the data arrays, find a place where another row may be stored. 5 Return PETSC_ERR_MEM if there is insufficient space to store the 6 row, so garbage collection and/or re-allocation may be done. 7 */ 8 #undef __FUNCT__ 9 #define __FUNCT__ "spbas_cholesky_row_alloc" 10 PetscErrorCode spbas_cholesky_row_alloc(spbas_matrix retval, PetscInt k, PetscInt r_nnz,PetscInt * n_alloc_used) 11 { 12 PetscFunctionBegin; 13 retval.icols[k] = &retval.alloc_icol[*n_alloc_used]; 14 retval.values[k] = &retval.alloc_val[*n_alloc_used]; 15 *n_alloc_used += r_nnz; 16 if (*n_alloc_used > retval.n_alloc_icol) 17 { 18 PetscFunctionReturn(PETSC_ERR_MEM); 19 } 20 else 21 { 22 PetscFunctionReturn(0); 23 } 24 } 25 26 27 /* 28 spbas_cholesky_garbage_collect: 29 move the rows which have been calculated so far, as well as 30 those currently under construction, to the front of the 31 array, while putting them in the proper order. 32 When it seems necessary, enlarge the current arrays. 33 34 NB: re-allocation is being simulated using 35 PetscMalloc, memcpy, PetscFree, because 36 PetscRealloc does not seem to exist. 37 38 */ 39 #undef __FUNCT__ 40 #define __FUNCT__ "spbas_cholesky_garbage_collect" 41 PetscErrorCode spbas_cholesky_garbage_collect( 42 spbas_matrix * result, // I/O: the Cholesky factor matrix being constructed. 43 // Only the storage, not the contents of this matrix 44 // is changed in this function 45 PetscInt i_row, // I : Number of rows for which the final contents are 46 // known 47 PetscInt *n_row_alloc_ok, // I/O: Number of rows which are already in their final 48 // places in the arrays: they need not be moved any more 49 PetscInt *n_alloc_used, // I/O: 50 PetscInt *max_row_nnz // I : Over-estimate of the number of nonzeros needed to 51 // store each row 52 ) 53 { 54 55 #undef garbage_debug 56 57 // PSEUDO-CODE: 58 // 1. Choose the appropriate size for the arrays 59 // 2. Rescue the arrays which would be lost during garbage collection 60 // 3. Reallocate and correct administration 61 // 4. Move all arrays so that they are in proper order 62 63 PetscInt i,j; 64 PetscInt nrows = result->nrows; 65 PetscInt n_alloc_ok=0; 66 PetscInt n_alloc_ok_max=0; 67 68 PetscErrorCode ierr; 69 PetscInt need_already= 0; 70 PetscInt n_rows_ahead=0; 71 PetscInt max_need_extra= 0; 72 PetscInt n_alloc_max, n_alloc_est, n_alloc; 73 PetscInt n_alloc_now = result->n_alloc_icol; 74 PetscInt *alloc_icol_old = result->alloc_icol; 75 PetscScalar *alloc_val_old = result->alloc_val; 76 PetscInt *icol_rescue; 77 PetscScalar *val_rescue; 78 PetscInt n_rescue; 79 PetscInt n_row_rescue; 80 PetscInt i_here, i_last, n_copy; 81 const PetscReal xtra_perc = 20; 82 83 PetscFunctionBegin; 84 85 //********************************************************* 86 // 1. Choose appropriate array size 87 // Count number of rows and memory usage which is already final 88 for (i=0;i<i_row; i++) 89 { 90 n_alloc_ok += result->row_nnz[i]; 91 n_alloc_ok_max += max_row_nnz[i]; 92 } 93 94 // Count currently needed memory usage and future memory requirements 95 // (max, predicted) 96 for (i=i_row; i<nrows; i++) 97 { 98 if (result->row_nnz[i]==0) 99 { 100 max_need_extra += max_row_nnz[i]; 101 } 102 else 103 { 104 need_already += max_row_nnz[i]; 105 n_rows_ahead++; 106 } 107 } 108 109 // Make maximal and realistic memory requirement estimates 110 n_alloc_max = n_alloc_ok + need_already + max_need_extra; 111 n_alloc_est = n_alloc_ok + need_already + (int) (((PetscReal) max_need_extra) * 112 ((PetscReal) n_alloc_ok) /((PetscReal) n_alloc_ok_max)); 113 114 // Choose array sizes 115 if (n_alloc_max == n_alloc_est) 116 { n_alloc = n_alloc_max; } 117 else if (n_alloc_now >= n_alloc_est) 118 { n_alloc = n_alloc_now; } 119 else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0)) 120 { n_alloc = n_alloc_max; } 121 else 122 { n_alloc = (int) (n_alloc_est * (1+xtra_perc/100.0)); } 123 124 // If new estimate is less than what we already have, 125 // don't reallocate, just garbage-collect 126 if (n_alloc_max != n_alloc_est && n_alloc < result->n_alloc_icol) 127 { 128 n_alloc = result->n_alloc_icol; 129 } 130 131 #ifdef garbage_debug 132 // Print allocation considerations 133 printf(" + Final version of rows 1..%d is available\n",i_row); 134 printf(" They have %d nonzeros: %6.2f per row\n",n_alloc_ok, 135 (PetscScalar) n_alloc_ok / (PetscScalar)i_row); 136 printf(" Without droptol, I would have had %d nonzeros\n", 137 n_alloc_ok_max); 138 printf(" So I use %6.2f%% of the allowed sparseness\n", 139 100.0 * (PetscScalar) n_alloc_ok / 140 (PetscScalar) n_alloc_ok_max); 141 printf(" + I have %d rows for which I need full sparseness pattern\n", 142 n_rows_ahead); 143 printf(" They require %d nonzeros\n", 144 need_already); 145 printf(" + After that, I shall need %d more rows\n", 146 nrows-n_rows_ahead-i_row); 147 printf(" Maximally, I will need %d nonzeros there\n", 148 max_need_extra); 149 printf(" Realistically, I expect to need %6.2f nonzeros there\n", 150 (PetscScalar) max_need_extra * 151 (PetscScalar) n_alloc_ok / 152 (PetscScalar) n_alloc_ok_max); 153 printf(" + New allocated arrays should be %d (max) or %d (realistic) long\n", 154 n_alloc_max, n_alloc_est); 155 156 #endif 157 // Motivate dimension choice 158 printf(" Allocating %d nonzeros: ",n_alloc); 159 if (n_alloc_max == n_alloc_est) 160 { printf("this is the correct size\n"); } 161 else if (n_alloc_now >= n_alloc_est) 162 { printf("the current size, which seems enough\n"); } 163 else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0)) 164 { printf("the maximum estimate\n"); } 165 else 166 { printf("%6.2f %% more than the estimate\n",xtra_perc); } 167 168 169 //********************************************************* 170 // 2. Rescue arrays which would be lost 171 // Count how many rows and nonzeros will have to be rescued 172 // when moving all arrays in place 173 n_row_rescue = 0; n_rescue = 0; 174 if (*n_row_alloc_ok==0) { *n_alloc_used = 0; } 175 else 176 { 177 i = *n_row_alloc_ok - 1; 178 *n_alloc_used = (result->icols[i]-result->alloc_icol) + 179 result->row_nnz[i]; 180 } 181 182 #ifdef garbage_debug 183 printf("\n\n"); 184 printf("Rows 0..%d are already OK\n", *n_row_alloc_ok); 185 printf("Rows 0..%d have the right dimensions\n", i_row-1); 186 #endif 187 188 for (i=*n_row_alloc_ok; i<nrows; i++) 189 { 190 i_here = result->icols[i]-result->alloc_icol; 191 i_last = i_here + result->row_nnz[i]; 192 if (result->row_nnz[i]>0) 193 { 194 if (*n_alloc_used > i_here || i_last > n_alloc) 195 { 196 n_rescue += result->row_nnz[i]; 197 n_row_rescue++; 198 //printf(" Rescue row %d: %d nonzeros\n",i,result->row_nnz[i]); 199 //printf(" (New start at %d, old start at %d)\n", 200 // *n_alloc_used, i_here); 201 } 202 203 if (i<i_row) { *n_alloc_used += result->row_nnz[i];} 204 else { *n_alloc_used += max_row_nnz[i];} 205 } 206 } 207 208 #ifdef garbage_debug 209 printf(" + %d arrays will have to be rescued: %d nnz\n", 210 n_row_rescue, n_rescue); 211 #endif 212 213 // Allocate rescue arrays 214 ierr = PetscMalloc( n_rescue * sizeof(int), &icol_rescue); 215 CHKERRQ(ierr); 216 ierr = PetscMalloc( n_rescue * sizeof(PetscScalar), &val_rescue); 217 CHKERRQ(ierr); 218 219 // Rescue the arrays which need rescuing 220 n_row_rescue = 0; n_rescue = 0; 221 if (*n_row_alloc_ok==0) { *n_alloc_used = 0; } 222 else 223 { 224 i = *n_row_alloc_ok - 1; 225 *n_alloc_used = (result->icols[i]-result->alloc_icol) + 226 result->row_nnz[i]; 227 } 228 229 for (i=*n_row_alloc_ok; i<nrows; i++) 230 { 231 i_here = result->icols[i]-result->alloc_icol; 232 i_last = i_here + result->row_nnz[i]; 233 if (result->row_nnz[i]>0) 234 { 235 if (*n_alloc_used > i_here || i_last > n_alloc) 236 { 237 //printf(" Rescuing row %d to array[%d..]\n",i,n_rescue); 238 memcpy((void*) &icol_rescue[n_rescue], 239 (void*) result->icols[i], 240 result->row_nnz[i] * sizeof(int)); 241 memcpy((void*) &val_rescue[n_rescue], 242 (void*) result->values[i], 243 result->row_nnz[i] * sizeof(PetscScalar)); 244 n_rescue += result->row_nnz[i]; 245 n_row_rescue++; 246 } 247 248 if (i<i_row) { *n_alloc_used += result->row_nnz[i];} 249 else { *n_alloc_used += max_row_nnz[i];} 250 } 251 } 252 253 //********************************************************* 254 // 3. Reallocate and correct administration 255 256 if (n_alloc != result->n_alloc_icol) 257 { 258 n_copy = MIN(n_alloc,result->n_alloc_icol); 259 260 // PETSC knows no REALLOC, so we'll REALLOC ourselves. 261 262 // Allocate new icol-data, copy old contents 263 ierr = PetscMalloc( n_alloc * sizeof(int), 264 &result->alloc_icol); CHKERRQ(ierr); 265 memcpy(result->alloc_icol, alloc_icol_old, n_copy*sizeof(int)); 266 267 // Update administration, Reset pointers to new arrays 268 result->n_alloc_icol = n_alloc; 269 for (i=0; i<nrows; i++) 270 { 271 result->icols[i] = result->alloc_icol + 272 (result->icols[i] - alloc_icol_old); 273 result->values[i] = result->alloc_val + 274 (result->icols[i] - result->alloc_icol); 275 } 276 277 // Delete old array 278 ierr = PetscFree(alloc_icol_old); CHKERRQ(ierr); 279 280 // Allocate new value-data, copy old contents 281 ierr = PetscMalloc( n_alloc * sizeof(PetscScalar), 282 &result->alloc_val); CHKERRQ(ierr); 283 memcpy(result->alloc_val, alloc_val_old, n_copy*sizeof(PetscScalar)); 284 285 // Update administration, Reset pointers to new arrays 286 result->n_alloc_val = n_alloc; 287 for (i=0; i<nrows; i++) 288 { 289 result->values[i] = result->alloc_val + 290 (result->icols[i] - result->alloc_icol); 291 } 292 293 // Delete old array 294 ierr = PetscFree(alloc_val_old); CHKERRQ(ierr); 295 } 296 297 298 //********************************************************* 299 // 4. Copy all the arrays to their proper places 300 n_row_rescue = 0; n_rescue = 0; 301 if (*n_row_alloc_ok==0) { *n_alloc_used = 0; } 302 else 303 { 304 i = *n_row_alloc_ok - 1; 305 *n_alloc_used = (result->icols[i]-result->alloc_icol) + 306 result->row_nnz[i]; 307 } 308 309 for (i=*n_row_alloc_ok; i<nrows; i++) 310 { 311 i_here = result->icols[i]-result->alloc_icol; 312 i_last = i_here + result->row_nnz[i]; 313 314 result->icols[i] = result->alloc_icol + *n_alloc_used; 315 result->values[i]= result->alloc_val + *n_alloc_used; 316 317 if (result->row_nnz[i]>0) 318 { 319 if (*n_alloc_used > i_here || i_last > n_alloc) 320 { 321 //printf(" Copying rescued array[%d..] to row %d\n", 322 // n_rescue,i); 323 //fflush(stdout); 324 memcpy((void*) result->icols[i], 325 (void*) &icol_rescue[n_rescue], 326 result->row_nnz[i] * sizeof(int)); 327 memcpy((void*) result->values[i], 328 (void*) &val_rescue[n_rescue], 329 result->row_nnz[i] * sizeof(PetscScalar)); 330 n_rescue += result->row_nnz[i]; 331 n_row_rescue++; 332 } 333 else 334 { 335 for (j=0; j<result->row_nnz[i]; j++) 336 { 337 result->icols[i][j] = result->alloc_icol[i_here+j]; 338 result->values[i][j] = result->alloc_val[ i_here+j]; 339 } 340 } 341 if (i<i_row) { *n_alloc_used += result->row_nnz[i];} 342 else { *n_alloc_used += max_row_nnz[i];} 343 } 344 } 345 346 // Delete the rescue arrays 347 ierr = PetscFree(icol_rescue); CHKERRQ(ierr); 348 ierr = PetscFree(val_rescue); CHKERRQ(ierr); 349 *n_row_alloc_ok = i_row; 350 351 PetscFunctionReturn(0); 352 353 } 354 355 /* 356 spbas_incomplete_cholesky: 357 incomplete Cholesky decomposition of a square, symmetric, 358 positive definite matrix. 359 360 In case negative diagonals are encountered, function returns 361 NEGATIVE_DIAGONAL. When this happens, call this function again 362 with a larger epsdiag_in, a less sparse pattern, and/or a smaller 363 droptol 364 */ 365 #undef __FUNCT__ 366 #define __FUNCT__ "spbas_incomplete_cholesky" 367 PetscErrorCode spbas_incomplete_cholesky( 368 Mat A, const PetscInt *rip, const PetscInt *riip, 369 spbas_matrix pattern, 370 PetscReal droptol, PetscReal epsdiag_in, 371 spbas_matrix * matrix_L) 372 { 373 374 #undef cholesky_debug 375 376 #ifdef cholesky_debug 377 PetscInt idebug=0; 378 const char * fmt= " (%d,%d) -> %12.8e\n"; 379 PetscScalar tmp; 380 #endif 381 382 PetscInt jL; 383 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 384 PetscInt *ai=a->i,*aj=a->j; 385 MatScalar *aa=a->a; 386 PetscInt nrows, ncols; 387 PetscInt *max_row_nnz; 388 PetscErrorCode ierr; 389 spbas_matrix retval; 390 PetscScalar * diag; 391 PetscScalar * val; 392 PetscScalar * lvec; 393 PetscScalar epsdiag; 394 PetscInt i,j,k; 395 const PetscTruth do_values = PETSC_TRUE; 396 PetscInt * r1_icol; PetscScalar *r1_val; 397 PetscInt * r_icol; PetscInt r_nnz; PetscScalar *r_val; 398 PetscInt * A_icol; PetscInt A_nnz; PetscScalar *A_val; 399 PetscInt * p_icol; PetscInt p_nnz; 400 PetscInt n_row_alloc_ok = 0; // number of rows which have been stored 401 // correctly in the matrix 402 PetscInt n_alloc_used = 0; // part of result->icols and result->values 403 // which is currently being used 404 PetscFunctionBegin; 405 406 // Convert the Manteuffel shift from 'fraction of average diagonal' to 407 // dimensioned value 408 ierr = MatGetSize(A, &nrows, &ncols); CHKERRQ(ierr); 409 ierr = MatGetTrace(A, &epsdiag); CHKERRQ(ierr); 410 epsdiag *= epsdiag_in / nrows; 411 ierr = PetscInfo2(PETSC_NULL," Dimensioned Manteuffel shift %G Drop tolerance %G\n", PetscRealPart(epsdiag),droptol);CHKERRQ(ierr); 412 413 if (droptol<1e-10) {droptol=1e-10;} 414 415 // Check dimensions 416 if ( (nrows != pattern.nrows) || 417 (ncols != pattern.ncols) || 418 (ncols != nrows) ) 419 { 420 SETERRQ( PETSC_ERR_ARG_INCOMP, 421 "Dimension error in spbas_incomplete_cholesky\n"); 422 } 423 424 // Set overall values 425 retval.nrows = nrows; 426 retval.ncols = nrows; 427 retval.nnz = pattern.nnz/10; 428 retval.col_idx_type = SPBAS_COLUMN_NUMBERS; 429 retval.block_data = PETSC_TRUE; 430 431 // Allocate sparseness pattern 432 ierr = spbas_allocate_pattern(&retval, do_values); CHKERRQ(ierr); 433 // Initial allocation of data 434 memset((void*) retval.row_nnz, 0, nrows*sizeof(int)); 435 ierr = spbas_allocate_data(&retval); CHKERRQ(ierr); 436 retval.nnz = 0; 437 438 // Allocate work arrays 439 ierr = PetscMalloc(nrows*sizeof(PetscScalar), &diag); CHKERRQ(ierr); 440 ierr = PetscMalloc(nrows*sizeof(PetscScalar), &val); CHKERRQ(ierr); 441 ierr = PetscMalloc(nrows*sizeof(PetscScalar), &lvec); CHKERRQ(ierr); 442 ierr = PetscMalloc(nrows*sizeof(int), &max_row_nnz); CHKERRQ(ierr); 443 if (!(diag && val && lvec && max_row_nnz)) 444 { 445 SETERRQ( PETSC_ERR_MEM, 446 "Allocation error in spbas_incomplete_cholesky\n"); 447 } 448 449 memset((void*) val, 0, nrows*sizeof(PetscScalar)); 450 memset((void*) lvec, 0, nrows*sizeof(PetscScalar)); 451 memset((void*) max_row_nnz, 0, nrows*sizeof(int)); 452 453 // Check correct format of sparseness pattern 454 if (pattern.col_idx_type != SPBAS_DIAGONAL_OFFSETS) 455 { 456 SETERRQ( PETSC_ERR_SUP_SYS, 457 "Error in spbas_incomplete_cholesky: must have diagonal offsets in pattern\n"); 458 } 459 460 // Count the nonzeros on transpose of pattern 461 for (i = 0; i<nrows; i++) 462 { 463 p_nnz = pattern.row_nnz[i]; 464 p_icol = pattern.icols[i]; 465 for (j=0; j<p_nnz; j++) 466 { 467 max_row_nnz[i+p_icol[j]]++; 468 } 469 } 470 471 // Calculate rows of L 472 for (i = 0; i<nrows; i++) 473 { 474 p_nnz = pattern.row_nnz[i]; 475 p_icol = pattern.icols[i]; 476 477 r_nnz = retval.row_nnz[i]; 478 r_icol = retval.icols[i]; 479 r_val = retval.values[i]; 480 A_nnz = ai[rip[i]+1] - ai[rip[i]]; 481 A_icol = &aj[ai[rip[i]]]; 482 A_val = &aa[ai[rip[i]]]; 483 484 // Calculate val += A(i,i:n)'; 485 for (j=0; j<A_nnz; j++) 486 { 487 k = riip[A_icol[j]]; 488 if (k>=i) { val[k] = A_val[j]; } 489 } 490 491 // Add regularization 492 val[i] += epsdiag; 493 494 // Calculate lvec = diag(D(0:i-1)) * L(0:i-1,i); 495 // val(i) = A(i,i) - L(0:i-1,i)' * lvec 496 for ( j=0; j<r_nnz; j++) 497 { 498 k = r_icol[j]; 499 lvec[k] = diag[k] * r_val[j]; 500 val[i] -= r_val[j] * lvec[k]; 501 } 502 503 // Calculate the new diagonal 504 diag[i] = val[i]; 505 if (PetscRealPart(diag[i])<droptol) 506 { 507 printf("Error in spbas_incomplete_cholesky:\n"); 508 printf("Negative diagonal in row %d\n",i+1); 509 510 // Delete the whole matrix at once. 511 ierr = spbas_delete(retval); 512 return NEGATIVE_DIAGONAL; 513 } 514 515 #ifdef cholesky_debug 516 if (idebug>=3) 517 { 518 printf("diagD * L (%d,:)=\n",i+1); 519 for (j=0; j<r_nnz; j++) 520 { 521 printf(fmt,r_icol[j]+1,i+1,lvec[r_icol[j]]); 522 } 523 printf("\n\n"); 524 } 525 526 if (idebug>=1) { printf("Diag = %12.8e\n\n\n",diag[i]);} 527 528 if (idebug>=2) 529 { 530 printf("L^T * diagD * L=\n"); 531 tmp = 0.0; 532 for ( j=0; j<r_nnz; j++) 533 { 534 k = r_icol[j]; 535 tmp += r_val[j] * lvec[k]; 536 } 537 printf(fmt, i+1,i+1,tmp); 538 } 539 #endif 540 541 // If necessary, allocate arrays 542 if (r_nnz==0) 543 { 544 ierr = spbas_cholesky_row_alloc( retval, i, 1, &n_alloc_used); 545 if (ierr == PETSC_ERR_MEM) 546 { 547 ierr = spbas_cholesky_garbage_collect( &retval, i, &n_row_alloc_ok, 548 &n_alloc_used, max_row_nnz); 549 CHKERRQ(ierr); 550 ierr = spbas_cholesky_row_alloc( retval, i, 1, &n_alloc_used); 551 } 552 CHKERRQ(ierr); 553 r_icol = retval.icols[i]; 554 r_val = retval.values[i]; 555 } 556 557 // Now, fill in 558 r_icol[r_nnz] = i; 559 r_val [r_nnz] = 1.0; 560 r_nnz++; 561 retval.row_nnz[i]++; 562 563 retval.nnz += r_nnz; 564 565 // Calculate 566 // val(i+1:n) = (A(i,i+1:n)- L(0:i-1,i+1:n)' * lvec)/diag(i) 567 for (j=1; j<p_nnz; j++) 568 { 569 k = i+p_icol[j]; 570 r1_icol = retval.icols[k]; 571 r1_val = retval.values[k]; 572 #ifdef cholesky_debug 573 tmp = 0.0; 574 #endif 575 for (jL=0; jL<retval.row_nnz[k]; jL++) 576 { 577 val[k] -= r1_val[jL] * lvec[r1_icol[jL]]; 578 #ifdef cholesky_debug 579 tmp += r1_val[jL] * lvec[r1_icol[jL]]; 580 #endif 581 } 582 #ifdef cholesky_debug 583 if (idebug>=2 && tmp!=0) {printf(fmt, k+1,i+1,tmp);} 584 #endif 585 val[k] /= diag[i]; 586 587 if (PetscAbsScalar(val[k]) > droptol || PetscAbsScalar(val[k])< -droptol) 588 { 589 // If necessary, allocate arrays 590 if (retval.row_nnz[k]==0) 591 { 592 ierr = spbas_cholesky_row_alloc( retval, k, max_row_nnz[k], 593 &n_alloc_used); 594 if (ierr == PETSC_ERR_MEM) 595 { 596 ierr = spbas_cholesky_garbage_collect( 597 &retval, i, &n_row_alloc_ok, 598 &n_alloc_used, max_row_nnz); 599 CHKERRQ(ierr); 600 ierr = spbas_cholesky_row_alloc( retval, k, max_row_nnz[k], 601 &n_alloc_used); 602 r_icol = retval.icols[i]; 603 r_val = retval.values[i]; 604 } 605 CHKERRQ(ierr); 606 } 607 608 // Now, fill in 609 retval.icols[k][retval.row_nnz[k]] = i; 610 retval.values[k][retval.row_nnz[k]] = val[k]; 611 retval.row_nnz[k]++; 612 } 613 val[k] = 0; 614 } 615 #ifdef cholesky_debug 616 if (idebug>=2) {printf("\n\n");} 617 #endif 618 619 // Erase the values used in the work arrays 620 for (j=0; j<r_nnz; j++) { lvec[r_icol[j]] = 0; } 621 622 #ifdef cholesky_debug 623 if (idebug>=2) 624 { 625 printf("new L=\n"); 626 for (jL=i; jL<nrows; jL++) 627 { 628 j = retval.row_nnz[jL]-1; 629 if (j>=0) 630 { 631 k = retval.icols[jL][j]; 632 if (k==i) printf(fmt, jL+1,i+1,retval.values[jL][j]); 633 } 634 } 635 printf("\n\n"); 636 } 637 #endif 638 } 639 640 ierr=PetscFree(lvec); CHKERRQ(ierr); ierr=PetscFree(val); CHKERRQ(ierr); 641 642 ierr = spbas_cholesky_garbage_collect( &retval, nrows, &n_row_alloc_ok, 643 &n_alloc_used, max_row_nnz); 644 ierr=PetscFree(max_row_nnz); CHKERRQ(ierr); 645 646 // Place the inverse of the diagonals in the matrix 647 for (i=0; i<nrows; i++) 648 { 649 r_nnz = retval.row_nnz[i]; 650 retval.values[i][r_nnz-1] = 1.0 / diag[i]; 651 for (j=0; j<r_nnz-1; j++) 652 { 653 retval.values[i][j] *= -1; 654 } 655 } 656 ierr=PetscFree(diag); CHKERRQ(ierr); 657 658 659 // Return successfully 660 *matrix_L = retval; 661 662 #ifdef cholesky_debug 663 for (i = 0; i<nrows; i++) 664 { 665 printf("Row %d\n",i+1); 666 for (j=0; j<retval.row_nnz[i]; j++) 667 { 668 printf(fmt, i+1,retval.icols[i][j]+1,retval.values[i][j]); 669 } 670 } 671 #endif 672 PetscFunctionReturn(0); 673 } 674 675