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 PetscInt need_already = 0; 49b2f1ab58SBarry Smith PetscInt n_rows_ahead =0; 50b2f1ab58SBarry Smith PetscInt max_need_extra= 0; 51b2f1ab58SBarry Smith PetscInt n_alloc_max, n_alloc_est, n_alloc; 52b2f1ab58SBarry Smith PetscInt n_alloc_now = result->n_alloc_icol; 53b2f1ab58SBarry Smith PetscInt *alloc_icol_old = result->alloc_icol; 54b2f1ab58SBarry Smith PetscScalar *alloc_val_old = result->alloc_val; 55b2f1ab58SBarry Smith PetscInt *icol_rescue; 56b2f1ab58SBarry Smith PetscScalar *val_rescue; 57b2f1ab58SBarry Smith PetscInt n_rescue; 58b2f1ab58SBarry Smith PetscInt n_row_rescue; 59b2f1ab58SBarry Smith PetscInt i_here, i_last, n_copy; 60d36aa34eSBarry Smith const PetscReal xtra_perc = 20; 61b2f1ab58SBarry Smith 62b2f1ab58SBarry Smith PetscFunctionBegin; 63c328eeadSBarry Smith /********************************************************* 64c328eeadSBarry Smith 1. Choose appropriate array size 65c328eeadSBarry Smith Count number of rows and memory usage which is already final */ 669ccc4a9bSBarry Smith for (i=0; i<i_row; i++) { 67b2f1ab58SBarry Smith n_alloc_ok += result->row_nnz[i]; 68b2f1ab58SBarry Smith n_alloc_ok_max += max_row_nnz[i]; 69b2f1ab58SBarry Smith } 70b2f1ab58SBarry Smith 71c328eeadSBarry Smith /* Count currently needed memory usage and future memory requirements 72c328eeadSBarry Smith (max, predicted)*/ 739ccc4a9bSBarry Smith for (i=i_row; i<nrows; i++) { 749ccc4a9bSBarry Smith if (!result->row_nnz[i]) { 75b2f1ab58SBarry Smith max_need_extra += max_row_nnz[i]; 769ccc4a9bSBarry Smith } else { 77b2f1ab58SBarry Smith need_already += max_row_nnz[i]; 78b2f1ab58SBarry Smith n_rows_ahead++; 79b2f1ab58SBarry Smith } 80b2f1ab58SBarry Smith } 81b2f1ab58SBarry Smith 82c328eeadSBarry Smith /* Make maximal and realistic memory requirement estimates */ 83b2f1ab58SBarry Smith n_alloc_max = n_alloc_ok + need_already + max_need_extra; 849ccc4a9bSBarry Smith n_alloc_est = n_alloc_ok + need_already + (int) (((PetscReal) max_need_extra) * ((PetscReal) n_alloc_ok) /((PetscReal) n_alloc_ok_max)); 85b2f1ab58SBarry Smith 86c328eeadSBarry Smith /* Choose array sizes */ 872205254eSKarl Rupp if (n_alloc_max == n_alloc_est) n_alloc = n_alloc_max; 882205254eSKarl Rupp else if (n_alloc_now >= n_alloc_est) n_alloc = n_alloc_now; 892205254eSKarl Rupp else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0)) n_alloc = n_alloc_max; 902205254eSKarl Rupp else n_alloc = (int) (n_alloc_est * (1+xtra_perc/100.0)); 91b2f1ab58SBarry Smith 92c328eeadSBarry Smith /* If new estimate is less than what we already have, 93c328eeadSBarry Smith don't reallocate, just garbage-collect */ 949ccc4a9bSBarry Smith if (n_alloc_max != n_alloc_est && n_alloc < result->n_alloc_icol) { 95b2f1ab58SBarry Smith n_alloc = result->n_alloc_icol; 96b2f1ab58SBarry Smith } 97b2f1ab58SBarry Smith 98c328eeadSBarry Smith /* Motivate dimension choice */ 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(NULL," Allocating %" PetscInt_FMT " nonzeros: ",n_alloc)); 1002205254eSKarl Rupp if (n_alloc_max == n_alloc_est) { 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(NULL,"this is the correct size\n")); 1022205254eSKarl Rupp } else if (n_alloc_now >= n_alloc_est) { 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(NULL,"the current size, which seems enough\n")); 1042205254eSKarl Rupp } else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0)) { 1055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(NULL,"the maximum estimate\n")); 1062205254eSKarl Rupp } else { 1075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(NULL,"%6.2f %% more than the estimate\n",(double)xtra_perc)); 1082205254eSKarl Rupp } 109b2f1ab58SBarry Smith 110c328eeadSBarry Smith /********************************************************** 111c328eeadSBarry Smith 2. Rescue arrays which would be lost 112c328eeadSBarry Smith Count how many rows and nonzeros will have to be rescued 113c328eeadSBarry Smith when moving all arrays in place */ 114b2f1ab58SBarry Smith n_row_rescue = 0; n_rescue = 0; 1152205254eSKarl Rupp if (*n_row_alloc_ok==0) *n_alloc_used = 0; 1169ccc4a9bSBarry Smith else { 117b2f1ab58SBarry Smith i = *n_row_alloc_ok - 1; 1182205254eSKarl Rupp 1199ccc4a9bSBarry Smith *n_alloc_used = (result->icols[i]-result->alloc_icol) + result->row_nnz[i]; 120b2f1ab58SBarry Smith } 121b2f1ab58SBarry Smith 1229ccc4a9bSBarry Smith for (i=*n_row_alloc_ok; i<nrows; i++) { 123b2f1ab58SBarry Smith i_here = result->icols[i]-result->alloc_icol; 124b2f1ab58SBarry Smith i_last = i_here + result->row_nnz[i]; 1259ccc4a9bSBarry Smith if (result->row_nnz[i]>0) { 1269ccc4a9bSBarry Smith if (*n_alloc_used > i_here || i_last > n_alloc) { 127b2f1ab58SBarry Smith n_rescue += result->row_nnz[i]; 128b2f1ab58SBarry Smith n_row_rescue++; 129b2f1ab58SBarry Smith } 130b2f1ab58SBarry Smith 1312205254eSKarl Rupp if (i<i_row) *n_alloc_used += result->row_nnz[i]; 1322205254eSKarl Rupp else *n_alloc_used += max_row_nnz[i]; 133b2f1ab58SBarry Smith } 134b2f1ab58SBarry Smith } 135b2f1ab58SBarry Smith 136c328eeadSBarry Smith /* Allocate rescue arrays */ 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n_rescue, &icol_rescue)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n_rescue, &val_rescue)); 139b2f1ab58SBarry Smith 140c328eeadSBarry Smith /* Rescue the arrays which need rescuing */ 141b2f1ab58SBarry Smith n_row_rescue = 0; n_rescue = 0; 1422205254eSKarl Rupp if (*n_row_alloc_ok==0) *n_alloc_used = 0; 1439ccc4a9bSBarry Smith else { 144b2f1ab58SBarry Smith i = *n_row_alloc_ok - 1; 1459ccc4a9bSBarry Smith *n_alloc_used = (result->icols[i]-result->alloc_icol) + result->row_nnz[i]; 146b2f1ab58SBarry Smith } 147b2f1ab58SBarry Smith 1489ccc4a9bSBarry Smith for (i=*n_row_alloc_ok; i<nrows; i++) { 149b2f1ab58SBarry Smith i_here = result->icols[i]-result->alloc_icol; 150b2f1ab58SBarry Smith i_last = i_here + result->row_nnz[i]; 1519ccc4a9bSBarry Smith if (result->row_nnz[i]>0) { 1524e1ff37aSBarry Smith if (*n_alloc_used > i_here || i_last > n_alloc) { 1535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(&icol_rescue[n_rescue], result->icols[i], result->row_nnz[i])); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(&val_rescue[n_rescue], result->values[i], result->row_nnz[i])); 155b2f1ab58SBarry Smith n_rescue += result->row_nnz[i]; 156b2f1ab58SBarry Smith n_row_rescue++; 157b2f1ab58SBarry Smith } 158b2f1ab58SBarry Smith 1592205254eSKarl Rupp if (i<i_row) *n_alloc_used += result->row_nnz[i]; 1602205254eSKarl Rupp else *n_alloc_used += max_row_nnz[i]; 161b2f1ab58SBarry Smith } 162b2f1ab58SBarry Smith } 163b2f1ab58SBarry Smith 164c328eeadSBarry Smith /********************************************************** 165c328eeadSBarry Smith 3. Reallocate and correct administration */ 166b2f1ab58SBarry Smith 1679ccc4a9bSBarry Smith if (n_alloc != result->n_alloc_icol) { 168cc6fc633SBarry Smith n_copy = PetscMin(n_alloc,result->n_alloc_icol); 169b2f1ab58SBarry Smith 170c328eeadSBarry Smith /* PETSC knows no REALLOC, so we'll REALLOC ourselves. 171b2f1ab58SBarry Smith 172c328eeadSBarry Smith Allocate new icol-data, copy old contents */ 1735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n_alloc, &result->alloc_icol)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(result->alloc_icol, alloc_icol_old, n_copy)); 175b2f1ab58SBarry Smith 176c328eeadSBarry Smith /* Update administration, Reset pointers to new arrays */ 177b2f1ab58SBarry Smith result->n_alloc_icol = n_alloc; 1789ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 1799ccc4a9bSBarry Smith result->icols[i] = result->alloc_icol + (result->icols[i] - alloc_icol_old); 1809ccc4a9bSBarry Smith result->values[i] = result->alloc_val + (result->icols[i] - result->alloc_icol); 181b2f1ab58SBarry Smith } 182b2f1ab58SBarry Smith 183c328eeadSBarry Smith /* Delete old array */ 1845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(alloc_icol_old)); 185b2f1ab58SBarry Smith 186c328eeadSBarry Smith /* Allocate new value-data, copy old contents */ 1875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n_alloc, &result->alloc_val)); 1885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(result->alloc_val, alloc_val_old, n_copy)); 189b2f1ab58SBarry Smith 190c328eeadSBarry Smith /* Update administration, Reset pointers to new arrays */ 191b2f1ab58SBarry Smith result->n_alloc_val = n_alloc; 1929ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 1939ccc4a9bSBarry Smith result->values[i] = result->alloc_val + (result->icols[i] - result->alloc_icol); 194b2f1ab58SBarry Smith } 195b2f1ab58SBarry Smith 196c328eeadSBarry Smith /* Delete old array */ 1975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(alloc_val_old)); 198b2f1ab58SBarry Smith } 199b2f1ab58SBarry Smith 200c328eeadSBarry Smith /********************************************************* 201c328eeadSBarry Smith 4. Copy all the arrays to their proper places */ 202b2f1ab58SBarry Smith n_row_rescue = 0; n_rescue = 0; 2032205254eSKarl Rupp if (*n_row_alloc_ok==0) *n_alloc_used = 0; 2049ccc4a9bSBarry Smith else { 205b2f1ab58SBarry Smith i = *n_row_alloc_ok - 1; 2062205254eSKarl Rupp 2079ccc4a9bSBarry Smith *n_alloc_used = (result->icols[i]-result->alloc_icol) + result->row_nnz[i]; 208b2f1ab58SBarry Smith } 209b2f1ab58SBarry Smith 2109ccc4a9bSBarry Smith for (i=*n_row_alloc_ok; i<nrows; i++) { 211b2f1ab58SBarry Smith i_here = result->icols[i]-result->alloc_icol; 212b2f1ab58SBarry Smith i_last = i_here + result->row_nnz[i]; 213b2f1ab58SBarry Smith 214b2f1ab58SBarry Smith result->icols[i] = result->alloc_icol + *n_alloc_used; 215b2f1ab58SBarry Smith result->values[i]= result->alloc_val + *n_alloc_used; 216b2f1ab58SBarry Smith 2179ccc4a9bSBarry Smith if (result->row_nnz[i]>0) { 2189ccc4a9bSBarry Smith if (*n_alloc_used > i_here || i_last > n_alloc) { 2195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(result->icols[i], &icol_rescue[n_rescue], result->row_nnz[i])); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(result->values[i],&val_rescue[n_rescue],result->row_nnz[i])); 2212205254eSKarl Rupp 222b2f1ab58SBarry Smith n_rescue += result->row_nnz[i]; 223b2f1ab58SBarry Smith n_row_rescue++; 2249ccc4a9bSBarry Smith } else { 2259ccc4a9bSBarry Smith for (j=0; j<result->row_nnz[i]; j++) { 226b2f1ab58SBarry Smith result->icols[i][j] = result->alloc_icol[i_here+j]; 227b2f1ab58SBarry Smith result->values[i][j] = result->alloc_val[i_here+j]; 228b2f1ab58SBarry Smith } 229b2f1ab58SBarry Smith } 2302205254eSKarl Rupp if (i<i_row) *n_alloc_used += result->row_nnz[i]; 2312205254eSKarl Rupp else *n_alloc_used += max_row_nnz[i]; 232b2f1ab58SBarry Smith } 233b2f1ab58SBarry Smith } 234b2f1ab58SBarry Smith 235c328eeadSBarry Smith /* Delete the rescue arrays */ 2365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(icol_rescue)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(val_rescue)); 2382205254eSKarl Rupp 239b2f1ab58SBarry Smith *n_row_alloc_ok = i_row; 240b2f1ab58SBarry Smith PetscFunctionReturn(0); 241b2f1ab58SBarry Smith } 242b2f1ab58SBarry Smith 243b2f1ab58SBarry Smith /* 244b2f1ab58SBarry Smith spbas_incomplete_cholesky: 245b2f1ab58SBarry Smith incomplete Cholesky decomposition of a square, symmetric, 246b2f1ab58SBarry Smith positive definite matrix. 247b2f1ab58SBarry Smith 248b2f1ab58SBarry Smith In case negative diagonals are encountered, function returns 249b2f1ab58SBarry Smith NEGATIVE_DIAGONAL. When this happens, call this function again 250b2f1ab58SBarry Smith with a larger epsdiag_in, a less sparse pattern, and/or a smaller 251b2f1ab58SBarry Smith droptol 252b2f1ab58SBarry Smith */ 253cc442ecdSBarry 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) 254b2f1ab58SBarry Smith { 255b2f1ab58SBarry Smith PetscInt jL; 256b2f1ab58SBarry Smith Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data; 257b2f1ab58SBarry Smith PetscInt *ai=a->i,*aj=a->j; 258b2f1ab58SBarry Smith MatScalar *aa=a->a; 259b2f1ab58SBarry Smith PetscInt nrows, ncols; 260b2f1ab58SBarry Smith PetscInt *max_row_nnz; 261b2f1ab58SBarry Smith spbas_matrix retval; 262b2f1ab58SBarry Smith PetscScalar *diag; 263b2f1ab58SBarry Smith PetscScalar *val; 264b2f1ab58SBarry Smith PetscScalar *lvec; 265b2f1ab58SBarry Smith PetscScalar epsdiag; 266b2f1ab58SBarry Smith PetscInt i,j,k; 267ace3abfcSBarry Smith const PetscBool do_values = PETSC_TRUE; 2689ccc4a9bSBarry Smith PetscInt *r1_icol; 2699ccc4a9bSBarry Smith PetscScalar *r1_val; 2709ccc4a9bSBarry Smith PetscInt *r_icol; 2719ccc4a9bSBarry Smith PetscInt r_nnz; 2729ccc4a9bSBarry Smith PetscScalar *r_val; 2739ccc4a9bSBarry Smith PetscInt *A_icol; 2749ccc4a9bSBarry Smith PetscInt A_nnz; 2759ccc4a9bSBarry Smith PetscScalar *A_val; 2769ccc4a9bSBarry Smith PetscInt *p_icol; 2779ccc4a9bSBarry Smith PetscInt p_nnz; 2789ccc4a9bSBarry Smith PetscInt n_row_alloc_ok = 0; /* number of rows which have been stored correctly in the matrix */ 2799ccc4a9bSBarry Smith PetscInt n_alloc_used = 0; /* part of result->icols and result->values which is currently being used */ 280b2f1ab58SBarry Smith 2819ccc4a9bSBarry Smith PetscFunctionBegin; 2824e1ff37aSBarry Smith /* Convert the Manteuffel shift from 'fraction of average diagonal' to dimensioned value */ 2835f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A, &nrows, &ncols)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetTrace(A, &epsdiag)); 2852205254eSKarl Rupp 286b2f1ab58SBarry Smith epsdiag *= epsdiag_in / nrows; 2872205254eSKarl Rupp 2885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(NULL," Dimensioned Manteuffel shift %g Drop tolerance %g\n", (double)PetscRealPart(epsdiag),(double)droptol)); 289b2f1ab58SBarry Smith 2902205254eSKarl Rupp if (droptol<1e-10) droptol=1e-10; 291b2f1ab58SBarry Smith 292b2f1ab58SBarry Smith retval.nrows = nrows; 293b2f1ab58SBarry Smith retval.ncols = nrows; 294b2f1ab58SBarry Smith retval.nnz = pattern.nnz/10; 295b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_COLUMN_NUMBERS; 296b2f1ab58SBarry Smith retval.block_data = PETSC_TRUE; 297b2f1ab58SBarry Smith 2985f80ce2aSJacob Faibussowitsch CHKERRQ(spbas_allocate_pattern(&retval, do_values)); 2995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(retval.row_nnz, nrows)); 3005f80ce2aSJacob Faibussowitsch CHKERRQ(spbas_allocate_data(&retval)); 301b2f1ab58SBarry Smith retval.nnz = 0; 302b2f1ab58SBarry Smith 3035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows, &diag)); 3045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(nrows, &val)); 3055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(nrows, &lvec)); 3065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(nrows, &max_row_nnz)); 307b2f1ab58SBarry Smith 308c328eeadSBarry Smith /* Count the nonzeros on transpose of pattern */ 3094e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 310b2f1ab58SBarry Smith p_nnz = pattern.row_nnz[i]; 311b2f1ab58SBarry Smith p_icol = pattern.icols[i]; 3124e1ff37aSBarry Smith for (j=0; j<p_nnz; j++) { 313b2f1ab58SBarry Smith max_row_nnz[i+p_icol[j]]++; 314b2f1ab58SBarry Smith } 315b2f1ab58SBarry Smith } 316b2f1ab58SBarry Smith 317c328eeadSBarry Smith /* Calculate rows of L */ 3184e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 319b2f1ab58SBarry Smith p_nnz = pattern.row_nnz[i]; 320b2f1ab58SBarry Smith p_icol = pattern.icols[i]; 321b2f1ab58SBarry Smith 322b2f1ab58SBarry Smith r_nnz = retval.row_nnz[i]; 323b2f1ab58SBarry Smith r_icol = retval.icols[i]; 324b2f1ab58SBarry Smith r_val = retval.values[i]; 325b2f1ab58SBarry Smith A_nnz = ai[rip[i]+1] - ai[rip[i]]; 326b2f1ab58SBarry Smith A_icol = &aj[ai[rip[i]]]; 327b2f1ab58SBarry Smith A_val = &aa[ai[rip[i]]]; 328b2f1ab58SBarry Smith 329c328eeadSBarry Smith /* Calculate val += A(i,i:n)'; */ 3304e1ff37aSBarry Smith for (j=0; j<A_nnz; j++) { 331b2f1ab58SBarry Smith k = riip[A_icol[j]]; 3322205254eSKarl Rupp if (k>=i) val[k] = A_val[j]; 333b2f1ab58SBarry Smith } 334b2f1ab58SBarry Smith 335c328eeadSBarry Smith /* Add regularization */ 336b2f1ab58SBarry Smith val[i] += epsdiag; 337b2f1ab58SBarry Smith 338c328eeadSBarry Smith /* Calculate lvec = diag(D(0:i-1)) * L(0:i-1,i); 339c328eeadSBarry Smith val(i) = A(i,i) - L(0:i-1,i)' * lvec */ 3404e1ff37aSBarry Smith for (j=0; j<r_nnz; j++) { 341b2f1ab58SBarry Smith k = r_icol[j]; 342b2f1ab58SBarry Smith lvec[k] = diag[k] * r_val[j]; 343b2f1ab58SBarry Smith val[i] -= r_val[j] * lvec[k]; 344b2f1ab58SBarry Smith } 345b2f1ab58SBarry Smith 346c328eeadSBarry Smith /* Calculate the new diagonal */ 347b2f1ab58SBarry Smith diag[i] = val[i]; 3484e1ff37aSBarry Smith if (PetscRealPart(diag[i])<droptol) { 3495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(NULL,"Error in spbas_incomplete_cholesky:\n")); 3505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(NULL,"Negative diagonal in row %" PetscInt_FMT "\n",i+1)); 351b2f1ab58SBarry Smith 352c328eeadSBarry Smith /* Delete the whole matrix at once. */ 3535f80ce2aSJacob Faibussowitsch CHKERRQ(spbas_delete(retval)); 354cc442ecdSBarry Smith *success = PETSC_FALSE; 355cc442ecdSBarry Smith PetscFunctionReturn(0); 356b2f1ab58SBarry Smith } 357b2f1ab58SBarry Smith 358c328eeadSBarry Smith /* If necessary, allocate arrays */ 3594e1ff37aSBarry Smith if (r_nnz==0) { 360cc442ecdSBarry Smith PetscBool success = spbas_cholesky_row_alloc(retval, i, 1, &n_alloc_used); 361*28b400f6SJacob Faibussowitsch PetscCheck(success,PETSC_COMM_SELF,PETSC_ERR_MEM,"spbas_cholesky_row_alloc() failed"); 362b2f1ab58SBarry Smith r_icol = retval.icols[i]; 363b2f1ab58SBarry Smith r_val = retval.values[i]; 364b2f1ab58SBarry Smith } 365b2f1ab58SBarry Smith 366c328eeadSBarry Smith /* Now, fill in */ 367b2f1ab58SBarry Smith r_icol[r_nnz] = i; 368b2f1ab58SBarry Smith r_val [r_nnz] = 1.0; 369b2f1ab58SBarry Smith r_nnz++; 370b2f1ab58SBarry Smith retval.row_nnz[i]++; 371b2f1ab58SBarry Smith 372b2f1ab58SBarry Smith retval.nnz += r_nnz; 373b2f1ab58SBarry Smith 374c328eeadSBarry Smith /* Calculate 375c328eeadSBarry Smith val(i+1:n) = (A(i,i+1:n)- L(0:i-1,i+1:n)' * lvec)/diag(i) */ 3764e1ff37aSBarry Smith for (j=1; j<p_nnz; j++) { 377b2f1ab58SBarry Smith k = i+p_icol[j]; 378b2f1ab58SBarry Smith r1_icol = retval.icols[k]; 379b2f1ab58SBarry Smith r1_val = retval.values[k]; 3804e1ff37aSBarry Smith for (jL=0; jL<retval.row_nnz[k]; jL++) { 381b2f1ab58SBarry Smith val[k] -= r1_val[jL] * lvec[r1_icol[jL]]; 382b2f1ab58SBarry Smith } 383b2f1ab58SBarry Smith val[k] /= diag[i]; 384b2f1ab58SBarry Smith 3854e1ff37aSBarry Smith if (PetscAbsScalar(val[k]) > droptol || PetscAbsScalar(val[k])< -droptol) { 386c328eeadSBarry Smith /* If necessary, allocate arrays */ 387cc442ecdSBarry Smith if (!retval.row_nnz[k]) { 388cc442ecdSBarry Smith PetscBool flag,success = spbas_cholesky_row_alloc(retval, k, max_row_nnz[k], &n_alloc_used); 389cc442ecdSBarry Smith if (!success) { 3905f80ce2aSJacob Faibussowitsch CHKERRQ(spbas_cholesky_garbage_collect(&retval, i, &n_row_alloc_ok, &n_alloc_used, max_row_nnz)); 391cc442ecdSBarry Smith flag = spbas_cholesky_row_alloc(retval, k, max_row_nnz[k], &n_alloc_used); 392*28b400f6SJacob Faibussowitsch PetscCheck(flag,PETSC_COMM_SELF,PETSC_ERR_MEM,"Allocation in spbas_cholesky_row_alloc() failed"); 393b2f1ab58SBarry Smith r_icol = retval.icols[i]; 394cc442ecdSBarry Smith } 395b2f1ab58SBarry Smith } 396b2f1ab58SBarry Smith 397b2f1ab58SBarry Smith retval.icols[k][retval.row_nnz[k]] = i; 398b2f1ab58SBarry Smith retval.values[k][retval.row_nnz[k]] = val[k]; 399b2f1ab58SBarry Smith retval.row_nnz[k]++; 400b2f1ab58SBarry Smith } 401b2f1ab58SBarry Smith val[k] = 0; 402b2f1ab58SBarry Smith } 403b2f1ab58SBarry Smith 40471d55d18SBarry Smith /* Erase the values used in the work arrays */ 4052205254eSKarl Rupp for (j=0; j<r_nnz; j++) lvec[r_icol[j]] = 0; 406b2f1ab58SBarry Smith } 407b2f1ab58SBarry Smith 4085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(lvec)); 4095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(val)); 410b2f1ab58SBarry Smith 4115f80ce2aSJacob Faibussowitsch CHKERRQ(spbas_cholesky_garbage_collect(&retval, nrows, &n_row_alloc_ok, &n_alloc_used, max_row_nnz)); 4125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(max_row_nnz)); 413b2f1ab58SBarry Smith 414c328eeadSBarry Smith /* Place the inverse of the diagonals in the matrix */ 4159ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 416b2f1ab58SBarry Smith r_nnz = retval.row_nnz[i]; 4172205254eSKarl Rupp 418b2f1ab58SBarry Smith retval.values[i][r_nnz-1] = 1.0 / diag[i]; 4199ccc4a9bSBarry Smith for (j=0; j<r_nnz-1; j++) { 420b2f1ab58SBarry Smith retval.values[i][j] *= -1; 421b2f1ab58SBarry Smith } 422b2f1ab58SBarry Smith } 4235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(diag)); 424b2f1ab58SBarry Smith *matrix_L = retval; 425cc442ecdSBarry Smith *success = PETSC_TRUE; 426b2f1ab58SBarry Smith PetscFunctionReturn(0); 427b2f1ab58SBarry Smith } 428