1b2f1ab58SBarry Smith 2b2f1ab58SBarry Smith /* 3b2f1ab58SBarry Smith spbas_cholesky_row_alloc: 4b2f1ab58SBarry Smith in the data arrays, find a place where another row may be stored. 5b2f1ab58SBarry Smith Return PETSC_ERR_MEM if there is insufficient space to store the 6b2f1ab58SBarry Smith row, so garbage collection and/or re-allocation may be done. 7b2f1ab58SBarry Smith */ 8cc442ecdSBarry Smith PetscBool spbas_cholesky_row_alloc(spbas_matrix retval, PetscInt k, PetscInt r_nnz,PetscInt * n_alloc_used) 9b2f1ab58SBarry Smith { 10b2f1ab58SBarry Smith PetscFunctionBegin; 11b2f1ab58SBarry Smith retval.icols[k] = &retval.alloc_icol[*n_alloc_used]; 12b2f1ab58SBarry Smith retval.values[k] = &retval.alloc_val[*n_alloc_used]; 13b2f1ab58SBarry Smith *n_alloc_used += r_nnz; 14cc442ecdSBarry Smith if (*n_alloc_used > retval.n_alloc_icol) PetscFunctionReturn(PETSC_FALSE); 15cc442ecdSBarry Smith PetscFunctionReturn(PETSC_TRUE); 16b2f1ab58SBarry Smith } 17b2f1ab58SBarry Smith 18b2f1ab58SBarry Smith /* 19b2f1ab58SBarry Smith spbas_cholesky_garbage_collect: 20b2f1ab58SBarry Smith move the rows which have been calculated so far, as well as 21b2f1ab58SBarry Smith those currently under construction, to the front of the 22b2f1ab58SBarry Smith array, while putting them in the proper order. 23b2f1ab58SBarry Smith When it seems necessary, enlarge the current arrays. 24b2f1ab58SBarry Smith 25b2f1ab58SBarry Smith NB: re-allocation is being simulated using 26b2f1ab58SBarry Smith PetscMalloc, memcpy, PetscFree, because 27b2f1ab58SBarry Smith PetscRealloc does not seem to exist. 28b2f1ab58SBarry Smith 29b2f1ab58SBarry Smith */ 30be332245SKarl Rupp PetscErrorCode spbas_cholesky_garbage_collect(spbas_matrix *result, /* I/O: the Cholesky factor matrix being constructed. 312205254eSKarl Rupp Only the storage, not the contents of this matrix is changed in this function */ 32be332245SKarl Rupp PetscInt i_row, /* I : Number of rows for which the final contents are known */ 33c328eeadSBarry Smith PetscInt *n_row_alloc_ok, /* I/O: Number of rows which are already in their final 34c328eeadSBarry Smith places in the arrays: they need not be moved any more */ 35c328eeadSBarry Smith PetscInt *n_alloc_used, /* I/O: */ 36be332245SKarl Rupp PetscInt *max_row_nnz) /* I : Over-estimate of the number of nonzeros needed to store each row */ 37b2f1ab58SBarry Smith { 38c328eeadSBarry Smith /* PSEUDO-CODE: 39c328eeadSBarry Smith 1. Choose the appropriate size for the arrays 40c328eeadSBarry Smith 2. Rescue the arrays which would be lost during garbage collection 41c328eeadSBarry Smith 3. Reallocate and correct administration 42c328eeadSBarry Smith 4. Move all arrays so that they are in proper order */ 43b2f1ab58SBarry Smith 44b2f1ab58SBarry Smith PetscInt i,j; 45b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 46b2f1ab58SBarry Smith PetscInt n_alloc_ok =0; 47b2f1ab58SBarry Smith PetscInt n_alloc_ok_max=0; 48b2f1ab58SBarry Smith PetscErrorCode ierr; 49b2f1ab58SBarry Smith PetscInt need_already = 0; 50b2f1ab58SBarry Smith PetscInt n_rows_ahead =0; 51b2f1ab58SBarry Smith PetscInt max_need_extra= 0; 52b2f1ab58SBarry Smith PetscInt n_alloc_max, n_alloc_est, n_alloc; 53b2f1ab58SBarry Smith PetscInt n_alloc_now = result->n_alloc_icol; 54b2f1ab58SBarry Smith PetscInt *alloc_icol_old = result->alloc_icol; 55b2f1ab58SBarry Smith PetscScalar *alloc_val_old = result->alloc_val; 56b2f1ab58SBarry Smith PetscInt *icol_rescue; 57b2f1ab58SBarry Smith PetscScalar *val_rescue; 58b2f1ab58SBarry Smith PetscInt n_rescue; 59b2f1ab58SBarry Smith PetscInt n_row_rescue; 60b2f1ab58SBarry Smith PetscInt i_here, i_last, n_copy; 61d36aa34eSBarry Smith const PetscReal xtra_perc = 20; 62b2f1ab58SBarry Smith 63b2f1ab58SBarry Smith PetscFunctionBegin; 64c328eeadSBarry Smith /********************************************************* 65c328eeadSBarry Smith 1. Choose appropriate array size 66c328eeadSBarry Smith Count number of rows and memory usage which is already final */ 679ccc4a9bSBarry Smith for (i=0; i<i_row; i++) { 68b2f1ab58SBarry Smith n_alloc_ok += result->row_nnz[i]; 69b2f1ab58SBarry Smith n_alloc_ok_max += max_row_nnz[i]; 70b2f1ab58SBarry Smith } 71b2f1ab58SBarry Smith 72c328eeadSBarry Smith /* Count currently needed memory usage and future memory requirements 73c328eeadSBarry Smith (max, predicted)*/ 749ccc4a9bSBarry Smith for (i=i_row; i<nrows; i++) { 759ccc4a9bSBarry Smith if (!result->row_nnz[i]) { 76b2f1ab58SBarry Smith max_need_extra += max_row_nnz[i]; 779ccc4a9bSBarry Smith } else { 78b2f1ab58SBarry Smith need_already += max_row_nnz[i]; 79b2f1ab58SBarry Smith n_rows_ahead++; 80b2f1ab58SBarry Smith } 81b2f1ab58SBarry Smith } 82b2f1ab58SBarry Smith 83c328eeadSBarry Smith /* Make maximal and realistic memory requirement estimates */ 84b2f1ab58SBarry Smith n_alloc_max = n_alloc_ok + need_already + max_need_extra; 859ccc4a9bSBarry Smith n_alloc_est = n_alloc_ok + need_already + (int) (((PetscReal) max_need_extra) * ((PetscReal) n_alloc_ok) /((PetscReal) n_alloc_ok_max)); 86b2f1ab58SBarry Smith 87c328eeadSBarry Smith /* Choose array sizes */ 882205254eSKarl Rupp if (n_alloc_max == n_alloc_est) n_alloc = n_alloc_max; 892205254eSKarl Rupp else if (n_alloc_now >= n_alloc_est) n_alloc = n_alloc_now; 902205254eSKarl Rupp else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0)) n_alloc = n_alloc_max; 912205254eSKarl Rupp else n_alloc = (int) (n_alloc_est * (1+xtra_perc/100.0)); 92b2f1ab58SBarry Smith 93c328eeadSBarry Smith /* If new estimate is less than what we already have, 94c328eeadSBarry Smith don't reallocate, just garbage-collect */ 959ccc4a9bSBarry Smith if (n_alloc_max != n_alloc_est && n_alloc < result->n_alloc_icol) { 96b2f1ab58SBarry Smith n_alloc = result->n_alloc_icol; 97b2f1ab58SBarry Smith } 98b2f1ab58SBarry Smith 99c328eeadSBarry Smith /* Motivate dimension choice */ 100*7d3de750SJacob Faibussowitsch ierr = PetscInfo(NULL," Allocating %" PetscInt_FMT " nonzeros: ",n_alloc);CHKERRQ(ierr); 1012205254eSKarl Rupp if (n_alloc_max == n_alloc_est) { 1020298fd71SBarry Smith ierr = PetscInfo(NULL,"this is the correct size\n");CHKERRQ(ierr); 1032205254eSKarl Rupp } else if (n_alloc_now >= n_alloc_est) { 1040298fd71SBarry Smith ierr = PetscInfo(NULL,"the current size, which seems enough\n");CHKERRQ(ierr); 1052205254eSKarl Rupp } else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0)) { 1060298fd71SBarry Smith ierr = PetscInfo(NULL,"the maximum estimate\n");CHKERRQ(ierr); 1072205254eSKarl Rupp } else { 108*7d3de750SJacob Faibussowitsch ierr = PetscInfo(NULL,"%6.2f %% more than the estimate\n",(double)xtra_perc);CHKERRQ(ierr); 1092205254eSKarl Rupp } 110b2f1ab58SBarry Smith 111c328eeadSBarry Smith /********************************************************** 112c328eeadSBarry Smith 2. Rescue arrays which would be lost 113c328eeadSBarry Smith Count how many rows and nonzeros will have to be rescued 114c328eeadSBarry Smith when moving all arrays in place */ 115b2f1ab58SBarry Smith n_row_rescue = 0; n_rescue = 0; 1162205254eSKarl Rupp if (*n_row_alloc_ok==0) *n_alloc_used = 0; 1179ccc4a9bSBarry Smith else { 118b2f1ab58SBarry Smith i = *n_row_alloc_ok - 1; 1192205254eSKarl Rupp 1209ccc4a9bSBarry Smith *n_alloc_used = (result->icols[i]-result->alloc_icol) + result->row_nnz[i]; 121b2f1ab58SBarry Smith } 122b2f1ab58SBarry Smith 1239ccc4a9bSBarry Smith for (i=*n_row_alloc_ok; i<nrows; i++) { 124b2f1ab58SBarry Smith i_here = result->icols[i]-result->alloc_icol; 125b2f1ab58SBarry Smith i_last = i_here + result->row_nnz[i]; 1269ccc4a9bSBarry Smith if (result->row_nnz[i]>0) { 1279ccc4a9bSBarry Smith if (*n_alloc_used > i_here || i_last > n_alloc) { 128b2f1ab58SBarry Smith n_rescue += result->row_nnz[i]; 129b2f1ab58SBarry Smith n_row_rescue++; 130b2f1ab58SBarry Smith } 131b2f1ab58SBarry Smith 1322205254eSKarl Rupp if (i<i_row) *n_alloc_used += result->row_nnz[i]; 1332205254eSKarl Rupp else *n_alloc_used += max_row_nnz[i]; 134b2f1ab58SBarry Smith } 135b2f1ab58SBarry Smith } 136b2f1ab58SBarry Smith 137c328eeadSBarry Smith /* Allocate rescue arrays */ 138785e854fSJed Brown ierr = PetscMalloc1(n_rescue, &icol_rescue);CHKERRQ(ierr); 139785e854fSJed Brown ierr = PetscMalloc1(n_rescue, &val_rescue);CHKERRQ(ierr); 140b2f1ab58SBarry Smith 141c328eeadSBarry Smith /* Rescue the arrays which need rescuing */ 142b2f1ab58SBarry Smith n_row_rescue = 0; n_rescue = 0; 1432205254eSKarl Rupp if (*n_row_alloc_ok==0) *n_alloc_used = 0; 1449ccc4a9bSBarry Smith else { 145b2f1ab58SBarry Smith i = *n_row_alloc_ok - 1; 1469ccc4a9bSBarry Smith *n_alloc_used = (result->icols[i]-result->alloc_icol) + result->row_nnz[i]; 147b2f1ab58SBarry Smith } 148b2f1ab58SBarry Smith 1499ccc4a9bSBarry Smith for (i=*n_row_alloc_ok; i<nrows; i++) { 150b2f1ab58SBarry Smith i_here = result->icols[i]-result->alloc_icol; 151b2f1ab58SBarry Smith i_last = i_here + result->row_nnz[i]; 1529ccc4a9bSBarry Smith if (result->row_nnz[i]>0) { 1534e1ff37aSBarry Smith if (*n_alloc_used > i_here || i_last > n_alloc) { 154580bdb30SBarry Smith ierr = PetscArraycpy(&icol_rescue[n_rescue], result->icols[i], result->row_nnz[i]);CHKERRQ(ierr); 155580bdb30SBarry Smith ierr = PetscArraycpy(&val_rescue[n_rescue], result->values[i], result->row_nnz[i]);CHKERRQ(ierr); 156b2f1ab58SBarry Smith n_rescue += result->row_nnz[i]; 157b2f1ab58SBarry Smith n_row_rescue++; 158b2f1ab58SBarry Smith } 159b2f1ab58SBarry Smith 1602205254eSKarl Rupp if (i<i_row) *n_alloc_used += result->row_nnz[i]; 1612205254eSKarl Rupp else *n_alloc_used += max_row_nnz[i]; 162b2f1ab58SBarry Smith } 163b2f1ab58SBarry Smith } 164b2f1ab58SBarry Smith 165c328eeadSBarry Smith /********************************************************** 166c328eeadSBarry Smith 3. Reallocate and correct administration */ 167b2f1ab58SBarry Smith 1689ccc4a9bSBarry Smith if (n_alloc != result->n_alloc_icol) { 169cc6fc633SBarry Smith n_copy = PetscMin(n_alloc,result->n_alloc_icol); 170b2f1ab58SBarry Smith 171c328eeadSBarry Smith /* PETSC knows no REALLOC, so we'll REALLOC ourselves. 172b2f1ab58SBarry Smith 173c328eeadSBarry Smith Allocate new icol-data, copy old contents */ 174785e854fSJed Brown ierr = PetscMalloc1(n_alloc, &result->alloc_icol);CHKERRQ(ierr); 175580bdb30SBarry Smith ierr = PetscArraycpy(result->alloc_icol, alloc_icol_old, n_copy);CHKERRQ(ierr); 176b2f1ab58SBarry Smith 177c328eeadSBarry Smith /* Update administration, Reset pointers to new arrays */ 178b2f1ab58SBarry Smith result->n_alloc_icol = n_alloc; 1799ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 1809ccc4a9bSBarry Smith result->icols[i] = result->alloc_icol + (result->icols[i] - alloc_icol_old); 1819ccc4a9bSBarry Smith result->values[i] = result->alloc_val + (result->icols[i] - result->alloc_icol); 182b2f1ab58SBarry Smith } 183b2f1ab58SBarry Smith 184c328eeadSBarry Smith /* Delete old array */ 185b2f1ab58SBarry Smith ierr = PetscFree(alloc_icol_old);CHKERRQ(ierr); 186b2f1ab58SBarry Smith 187c328eeadSBarry Smith /* Allocate new value-data, copy old contents */ 188785e854fSJed Brown ierr = PetscMalloc1(n_alloc, &result->alloc_val);CHKERRQ(ierr); 189580bdb30SBarry Smith ierr = PetscArraycpy(result->alloc_val, alloc_val_old, n_copy);CHKERRQ(ierr); 190b2f1ab58SBarry Smith 191c328eeadSBarry Smith /* Update administration, Reset pointers to new arrays */ 192b2f1ab58SBarry Smith result->n_alloc_val = n_alloc; 1939ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 1949ccc4a9bSBarry Smith result->values[i] = result->alloc_val + (result->icols[i] - result->alloc_icol); 195b2f1ab58SBarry Smith } 196b2f1ab58SBarry Smith 197c328eeadSBarry Smith /* Delete old array */ 198b2f1ab58SBarry Smith ierr = PetscFree(alloc_val_old);CHKERRQ(ierr); 199b2f1ab58SBarry Smith } 200b2f1ab58SBarry Smith 201c328eeadSBarry Smith /********************************************************* 202c328eeadSBarry Smith 4. Copy all the arrays to their proper places */ 203b2f1ab58SBarry Smith n_row_rescue = 0; n_rescue = 0; 2042205254eSKarl Rupp if (*n_row_alloc_ok==0) *n_alloc_used = 0; 2059ccc4a9bSBarry Smith else { 206b2f1ab58SBarry Smith i = *n_row_alloc_ok - 1; 2072205254eSKarl Rupp 2089ccc4a9bSBarry Smith *n_alloc_used = (result->icols[i]-result->alloc_icol) + result->row_nnz[i]; 209b2f1ab58SBarry Smith } 210b2f1ab58SBarry Smith 2119ccc4a9bSBarry Smith for (i=*n_row_alloc_ok; i<nrows; i++) { 212b2f1ab58SBarry Smith i_here = result->icols[i]-result->alloc_icol; 213b2f1ab58SBarry Smith i_last = i_here + result->row_nnz[i]; 214b2f1ab58SBarry Smith 215b2f1ab58SBarry Smith result->icols[i] = result->alloc_icol + *n_alloc_used; 216b2f1ab58SBarry Smith result->values[i]= result->alloc_val + *n_alloc_used; 217b2f1ab58SBarry Smith 2189ccc4a9bSBarry Smith if (result->row_nnz[i]>0) { 2199ccc4a9bSBarry Smith if (*n_alloc_used > i_here || i_last > n_alloc) { 220580bdb30SBarry Smith ierr = PetscArraycpy(result->icols[i], &icol_rescue[n_rescue], result->row_nnz[i]);CHKERRQ(ierr); 221580bdb30SBarry Smith ierr = PetscArraycpy(result->values[i],&val_rescue[n_rescue],result->row_nnz[i]);CHKERRQ(ierr); 2222205254eSKarl Rupp 223b2f1ab58SBarry Smith n_rescue += result->row_nnz[i]; 224b2f1ab58SBarry Smith n_row_rescue++; 2259ccc4a9bSBarry Smith } else { 2269ccc4a9bSBarry Smith for (j=0; j<result->row_nnz[i]; j++) { 227b2f1ab58SBarry Smith result->icols[i][j] = result->alloc_icol[i_here+j]; 228b2f1ab58SBarry Smith result->values[i][j] = result->alloc_val[i_here+j]; 229b2f1ab58SBarry Smith } 230b2f1ab58SBarry Smith } 2312205254eSKarl Rupp if (i<i_row) *n_alloc_used += result->row_nnz[i]; 2322205254eSKarl Rupp else *n_alloc_used += max_row_nnz[i]; 233b2f1ab58SBarry Smith } 234b2f1ab58SBarry Smith } 235b2f1ab58SBarry Smith 236c328eeadSBarry Smith /* Delete the rescue arrays */ 237b2f1ab58SBarry Smith ierr = PetscFree(icol_rescue);CHKERRQ(ierr); 238b2f1ab58SBarry Smith ierr = PetscFree(val_rescue);CHKERRQ(ierr); 2392205254eSKarl Rupp 240b2f1ab58SBarry Smith *n_row_alloc_ok = i_row; 241b2f1ab58SBarry Smith PetscFunctionReturn(0); 242b2f1ab58SBarry Smith } 243b2f1ab58SBarry Smith 244b2f1ab58SBarry Smith /* 245b2f1ab58SBarry Smith spbas_incomplete_cholesky: 246b2f1ab58SBarry Smith incomplete Cholesky decomposition of a square, symmetric, 247b2f1ab58SBarry Smith positive definite matrix. 248b2f1ab58SBarry Smith 249b2f1ab58SBarry Smith In case negative diagonals are encountered, function returns 250b2f1ab58SBarry Smith NEGATIVE_DIAGONAL. When this happens, call this function again 251b2f1ab58SBarry Smith with a larger epsdiag_in, a less sparse pattern, and/or a smaller 252b2f1ab58SBarry Smith droptol 253b2f1ab58SBarry Smith */ 254cc442ecdSBarry Smith PetscErrorCode spbas_incomplete_cholesky(Mat A, const PetscInt *rip, const PetscInt *riip, spbas_matrix pattern, PetscReal droptol, PetscReal epsdiag_in, spbas_matrix * matrix_L,PetscBool *success) 255b2f1ab58SBarry Smith { 256b2f1ab58SBarry Smith PetscInt jL; 257b2f1ab58SBarry Smith Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data; 258b2f1ab58SBarry Smith PetscInt *ai=a->i,*aj=a->j; 259b2f1ab58SBarry Smith MatScalar *aa=a->a; 260b2f1ab58SBarry Smith PetscInt nrows, ncols; 261b2f1ab58SBarry Smith PetscInt *max_row_nnz; 262b2f1ab58SBarry Smith PetscErrorCode ierr; 263b2f1ab58SBarry Smith spbas_matrix retval; 264b2f1ab58SBarry Smith PetscScalar *diag; 265b2f1ab58SBarry Smith PetscScalar *val; 266b2f1ab58SBarry Smith PetscScalar *lvec; 267b2f1ab58SBarry Smith PetscScalar epsdiag; 268b2f1ab58SBarry Smith PetscInt i,j,k; 269ace3abfcSBarry Smith const PetscBool do_values = PETSC_TRUE; 2709ccc4a9bSBarry Smith PetscInt *r1_icol; 2719ccc4a9bSBarry Smith PetscScalar *r1_val; 2729ccc4a9bSBarry Smith PetscInt *r_icol; 2739ccc4a9bSBarry Smith PetscInt r_nnz; 2749ccc4a9bSBarry Smith PetscScalar *r_val; 2759ccc4a9bSBarry Smith PetscInt *A_icol; 2769ccc4a9bSBarry Smith PetscInt A_nnz; 2779ccc4a9bSBarry Smith PetscScalar *A_val; 2789ccc4a9bSBarry Smith PetscInt *p_icol; 2799ccc4a9bSBarry Smith PetscInt p_nnz; 2809ccc4a9bSBarry Smith PetscInt n_row_alloc_ok = 0; /* number of rows which have been stored correctly in the matrix */ 2819ccc4a9bSBarry Smith PetscInt n_alloc_used = 0; /* part of result->icols and result->values which is currently being used */ 282b2f1ab58SBarry Smith 2839ccc4a9bSBarry Smith PetscFunctionBegin; 2844e1ff37aSBarry Smith /* Convert the Manteuffel shift from 'fraction of average diagonal' to dimensioned value */ 285b2f1ab58SBarry Smith ierr = MatGetSize(A, &nrows, &ncols);CHKERRQ(ierr); 286b2f1ab58SBarry Smith ierr = MatGetTrace(A, &epsdiag);CHKERRQ(ierr); 2872205254eSKarl Rupp 288b2f1ab58SBarry Smith epsdiag *= epsdiag_in / nrows; 2892205254eSKarl Rupp 290*7d3de750SJacob Faibussowitsch ierr = PetscInfo(NULL," Dimensioned Manteuffel shift %g Drop tolerance %g\n", (double)PetscRealPart(epsdiag),(double)droptol);CHKERRQ(ierr); 291b2f1ab58SBarry Smith 2922205254eSKarl Rupp if (droptol<1e-10) droptol=1e-10; 293b2f1ab58SBarry Smith 294b2f1ab58SBarry Smith retval.nrows = nrows; 295b2f1ab58SBarry Smith retval.ncols = nrows; 296b2f1ab58SBarry Smith retval.nnz = pattern.nnz/10; 297b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_COLUMN_NUMBERS; 298b2f1ab58SBarry Smith retval.block_data = PETSC_TRUE; 299b2f1ab58SBarry Smith 300b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(&retval, do_values);CHKERRQ(ierr); 301580bdb30SBarry Smith ierr = PetscArrayzero(retval.row_nnz, nrows);CHKERRQ(ierr); 302b2f1ab58SBarry Smith ierr = spbas_allocate_data(&retval);CHKERRQ(ierr); 303b2f1ab58SBarry Smith retval.nnz = 0; 304b2f1ab58SBarry Smith 305785e854fSJed Brown ierr = PetscMalloc1(nrows, &diag);CHKERRQ(ierr); 306580bdb30SBarry Smith ierr = PetscCalloc1(nrows, &val);CHKERRQ(ierr); 307580bdb30SBarry Smith ierr = PetscCalloc1(nrows, &lvec);CHKERRQ(ierr); 308580bdb30SBarry Smith ierr = PetscCalloc1(nrows, &max_row_nnz);CHKERRQ(ierr); 309b2f1ab58SBarry Smith 310c328eeadSBarry Smith /* Count the nonzeros on transpose of pattern */ 3114e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 312b2f1ab58SBarry Smith p_nnz = pattern.row_nnz[i]; 313b2f1ab58SBarry Smith p_icol = pattern.icols[i]; 3144e1ff37aSBarry Smith for (j=0; j<p_nnz; j++) { 315b2f1ab58SBarry Smith max_row_nnz[i+p_icol[j]]++; 316b2f1ab58SBarry Smith } 317b2f1ab58SBarry Smith } 318b2f1ab58SBarry Smith 319c328eeadSBarry Smith /* Calculate rows of L */ 3204e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 321b2f1ab58SBarry Smith p_nnz = pattern.row_nnz[i]; 322b2f1ab58SBarry Smith p_icol = pattern.icols[i]; 323b2f1ab58SBarry Smith 324b2f1ab58SBarry Smith r_nnz = retval.row_nnz[i]; 325b2f1ab58SBarry Smith r_icol = retval.icols[i]; 326b2f1ab58SBarry Smith r_val = retval.values[i]; 327b2f1ab58SBarry Smith A_nnz = ai[rip[i]+1] - ai[rip[i]]; 328b2f1ab58SBarry Smith A_icol = &aj[ai[rip[i]]]; 329b2f1ab58SBarry Smith A_val = &aa[ai[rip[i]]]; 330b2f1ab58SBarry Smith 331c328eeadSBarry Smith /* Calculate val += A(i,i:n)'; */ 3324e1ff37aSBarry Smith for (j=0; j<A_nnz; j++) { 333b2f1ab58SBarry Smith k = riip[A_icol[j]]; 3342205254eSKarl Rupp if (k>=i) val[k] = A_val[j]; 335b2f1ab58SBarry Smith } 336b2f1ab58SBarry Smith 337c328eeadSBarry Smith /* Add regularization */ 338b2f1ab58SBarry Smith val[i] += epsdiag; 339b2f1ab58SBarry Smith 340c328eeadSBarry Smith /* Calculate lvec = diag(D(0:i-1)) * L(0:i-1,i); 341c328eeadSBarry Smith val(i) = A(i,i) - L(0:i-1,i)' * lvec */ 3424e1ff37aSBarry Smith for (j=0; j<r_nnz; j++) { 343b2f1ab58SBarry Smith k = r_icol[j]; 344b2f1ab58SBarry Smith lvec[k] = diag[k] * r_val[j]; 345b2f1ab58SBarry Smith val[i] -= r_val[j] * lvec[k]; 346b2f1ab58SBarry Smith } 347b2f1ab58SBarry Smith 348c328eeadSBarry Smith /* Calculate the new diagonal */ 349b2f1ab58SBarry Smith diag[i] = val[i]; 3504e1ff37aSBarry Smith if (PetscRealPart(diag[i])<droptol) { 351302440fdSBarry Smith ierr = PetscInfo(NULL,"Error in spbas_incomplete_cholesky:\n");CHKERRQ(ierr); 352*7d3de750SJacob Faibussowitsch ierr = PetscInfo(NULL,"Negative diagonal in row %" PetscInt_FMT "\n",i+1);CHKERRQ(ierr); 353b2f1ab58SBarry Smith 354c328eeadSBarry Smith /* Delete the whole matrix at once. */ 355302440fdSBarry Smith ierr = spbas_delete(retval);CHKERRQ(ierr); 356cc442ecdSBarry Smith *success = PETSC_FALSE; 357cc442ecdSBarry Smith PetscFunctionReturn(0); 358b2f1ab58SBarry Smith } 359b2f1ab58SBarry Smith 360c328eeadSBarry Smith /* If necessary, allocate arrays */ 3614e1ff37aSBarry Smith if (r_nnz==0) { 362cc442ecdSBarry Smith PetscBool success = spbas_cholesky_row_alloc(retval, i, 1, &n_alloc_used); 363cc442ecdSBarry Smith if (!success) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"spbas_cholesky_row_alloc() failed"); 364b2f1ab58SBarry Smith r_icol = retval.icols[i]; 365b2f1ab58SBarry Smith r_val = retval.values[i]; 366b2f1ab58SBarry Smith } 367b2f1ab58SBarry Smith 368c328eeadSBarry Smith /* Now, fill in */ 369b2f1ab58SBarry Smith r_icol[r_nnz] = i; 370b2f1ab58SBarry Smith r_val [r_nnz] = 1.0; 371b2f1ab58SBarry Smith r_nnz++; 372b2f1ab58SBarry Smith retval.row_nnz[i]++; 373b2f1ab58SBarry Smith 374b2f1ab58SBarry Smith retval.nnz += r_nnz; 375b2f1ab58SBarry Smith 376c328eeadSBarry Smith /* Calculate 377c328eeadSBarry Smith val(i+1:n) = (A(i,i+1:n)- L(0:i-1,i+1:n)' * lvec)/diag(i) */ 3784e1ff37aSBarry Smith for (j=1; j<p_nnz; j++) { 379b2f1ab58SBarry Smith k = i+p_icol[j]; 380b2f1ab58SBarry Smith r1_icol = retval.icols[k]; 381b2f1ab58SBarry Smith r1_val = retval.values[k]; 3824e1ff37aSBarry Smith for (jL=0; jL<retval.row_nnz[k]; jL++) { 383b2f1ab58SBarry Smith val[k] -= r1_val[jL] * lvec[r1_icol[jL]]; 384b2f1ab58SBarry Smith } 385b2f1ab58SBarry Smith val[k] /= diag[i]; 386b2f1ab58SBarry Smith 3874e1ff37aSBarry Smith if (PetscAbsScalar(val[k]) > droptol || PetscAbsScalar(val[k])< -droptol) { 388c328eeadSBarry Smith /* If necessary, allocate arrays */ 389cc442ecdSBarry Smith if (!retval.row_nnz[k]) { 390cc442ecdSBarry Smith PetscBool flag,success = spbas_cholesky_row_alloc(retval, k, max_row_nnz[k], &n_alloc_used); 391cc442ecdSBarry Smith if (!success) { 3929ccc4a9bSBarry Smith ierr = spbas_cholesky_garbage_collect(&retval, i, &n_row_alloc_ok, &n_alloc_used, max_row_nnz);CHKERRQ(ierr); 393cc442ecdSBarry Smith flag = spbas_cholesky_row_alloc(retval, k, max_row_nnz[k], &n_alloc_used); 394cc442ecdSBarry Smith if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Allocation in spbas_cholesky_row_alloc() failed"); 395b2f1ab58SBarry Smith r_icol = retval.icols[i]; 396cc442ecdSBarry Smith } 397b2f1ab58SBarry Smith } 398b2f1ab58SBarry Smith 399b2f1ab58SBarry Smith retval.icols[k][retval.row_nnz[k]] = i; 400b2f1ab58SBarry Smith retval.values[k][retval.row_nnz[k]] = val[k]; 401b2f1ab58SBarry Smith retval.row_nnz[k]++; 402b2f1ab58SBarry Smith } 403b2f1ab58SBarry Smith val[k] = 0; 404b2f1ab58SBarry Smith } 405b2f1ab58SBarry Smith 40671d55d18SBarry Smith /* Erase the values used in the work arrays */ 4072205254eSKarl Rupp for (j=0; j<r_nnz; j++) lvec[r_icol[j]] = 0; 408b2f1ab58SBarry Smith } 409b2f1ab58SBarry Smith 4109ccc4a9bSBarry Smith ierr = PetscFree(lvec);CHKERRQ(ierr); 4119ccc4a9bSBarry Smith ierr = PetscFree(val);CHKERRQ(ierr); 412b2f1ab58SBarry Smith 4139ccc4a9bSBarry Smith ierr = spbas_cholesky_garbage_collect(&retval, nrows, &n_row_alloc_ok, &n_alloc_used, max_row_nnz);CHKERRQ(ierr); 414b2f1ab58SBarry Smith ierr = PetscFree(max_row_nnz);CHKERRQ(ierr); 415b2f1ab58SBarry Smith 416c328eeadSBarry Smith /* Place the inverse of the diagonals in the matrix */ 4179ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 418b2f1ab58SBarry Smith r_nnz = retval.row_nnz[i]; 4192205254eSKarl Rupp 420b2f1ab58SBarry Smith retval.values[i][r_nnz-1] = 1.0 / diag[i]; 4219ccc4a9bSBarry Smith for (j=0; j<r_nnz-1; j++) { 422b2f1ab58SBarry Smith retval.values[i][j] *= -1; 423b2f1ab58SBarry Smith } 424b2f1ab58SBarry Smith } 425b2f1ab58SBarry Smith ierr = PetscFree(diag);CHKERRQ(ierr); 426b2f1ab58SBarry Smith *matrix_L = retval; 427cc442ecdSBarry Smith *success = PETSC_TRUE; 428b2f1ab58SBarry Smith PetscFunctionReturn(0); 429b2f1ab58SBarry Smith } 430