xref: /petsc/src/mat/impls/aij/seq/bas/spbas_cholesky.h (revision 11cc89d2bae31c43c795890e5d9b25a12b75f4f2)
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 */
89371c9d4SSatish Balay PetscBool spbas_cholesky_row_alloc(spbas_matrix retval, PetscInt k, PetscInt r_nnz, PetscInt *n_alloc_used) {
9b2f1ab58SBarry Smith   retval.icols[k]  = &retval.alloc_icol[*n_alloc_used];
10b2f1ab58SBarry Smith   retval.values[k] = &retval.alloc_val[*n_alloc_used];
11b2f1ab58SBarry Smith   *n_alloc_used += r_nnz;
12*11cc89d2SBarry Smith   return (*n_alloc_used > retval.n_alloc_icol) ? PETSC_FALSE : PETSC_TRUE;
13b2f1ab58SBarry Smith }
14b2f1ab58SBarry Smith 
15b2f1ab58SBarry Smith /*
16b2f1ab58SBarry Smith   spbas_cholesky_garbage_collect:
17b2f1ab58SBarry Smith      move the rows which have been calculated so far, as well as
18b2f1ab58SBarry Smith      those currently under construction, to the front of the
19b2f1ab58SBarry Smith      array, while putting them in the proper order.
20b2f1ab58SBarry Smith      When it seems necessary, enlarge the current arrays.
21b2f1ab58SBarry Smith 
22b2f1ab58SBarry Smith      NB: re-allocation is being simulated using
23b2f1ab58SBarry Smith          PetscMalloc, memcpy, PetscFree, because
24b2f1ab58SBarry Smith          PetscRealloc does not seem to exist.
25b2f1ab58SBarry Smith 
26b2f1ab58SBarry Smith */
27be332245SKarl Rupp PetscErrorCode spbas_cholesky_garbage_collect(spbas_matrix *result,         /* I/O: the Cholesky factor matrix being constructed.
282205254eSKarl Rupp                                                                                     Only the storage, not the contents of this matrix is changed in this function */
29be332245SKarl Rupp                                               PetscInt      i_row,          /* I  : Number of rows for which the final contents are known */
30c328eeadSBarry Smith                                               PetscInt     *n_row_alloc_ok, /* I/O: Number of rows which are already in their final
31c328eeadSBarry Smith                                                                                     places in the arrays: they need not be moved any more */
32c328eeadSBarry Smith                                               PetscInt     *n_alloc_used,   /* I/O:  */
33be332245SKarl Rupp                                               PetscInt     *max_row_nnz)        /* I  : Over-estimate of the number of nonzeros needed to store each row */
34b2f1ab58SBarry Smith {
35c328eeadSBarry Smith   /* PSEUDO-CODE:
36c328eeadSBarry Smith   1. Choose the appropriate size for the arrays
37c328eeadSBarry Smith   2. Rescue the arrays which would be lost during garbage collection
38c328eeadSBarry Smith   3. Reallocate and correct administration
39c328eeadSBarry Smith   4. Move all arrays so that they are in proper order */
40b2f1ab58SBarry Smith 
41b2f1ab58SBarry Smith   PetscInt        i, j;
42b2f1ab58SBarry Smith   PetscInt        nrows          = result->nrows;
43b2f1ab58SBarry Smith   PetscInt        n_alloc_ok     = 0;
44b2f1ab58SBarry Smith   PetscInt        n_alloc_ok_max = 0;
45b2f1ab58SBarry Smith   PetscInt        need_already   = 0;
46b2f1ab58SBarry Smith   PetscInt        max_need_extra = 0;
47b2f1ab58SBarry Smith   PetscInt        n_alloc_max, n_alloc_est, n_alloc;
48b2f1ab58SBarry Smith   PetscInt        n_alloc_now    = result->n_alloc_icol;
49b2f1ab58SBarry Smith   PetscInt       *alloc_icol_old = result->alloc_icol;
50b2f1ab58SBarry Smith   PetscScalar    *alloc_val_old  = result->alloc_val;
51b2f1ab58SBarry Smith   PetscInt       *icol_rescue;
52b2f1ab58SBarry Smith   PetscScalar    *val_rescue;
53b2f1ab58SBarry Smith   PetscInt        n_rescue;
54b2f1ab58SBarry Smith   PetscInt        i_here, i_last, n_copy;
55d36aa34eSBarry Smith   const PetscReal xtra_perc = 20;
56b2f1ab58SBarry Smith 
57b2f1ab58SBarry Smith   PetscFunctionBegin;
58c328eeadSBarry Smith   /*********************************************************
59c328eeadSBarry Smith   1. Choose appropriate array size
60c328eeadSBarry Smith   Count number of rows and memory usage which is already final */
619ccc4a9bSBarry Smith   for (i = 0; i < i_row; i++) {
62b2f1ab58SBarry Smith     n_alloc_ok += result->row_nnz[i];
63b2f1ab58SBarry Smith     n_alloc_ok_max += max_row_nnz[i];
64b2f1ab58SBarry Smith   }
65b2f1ab58SBarry Smith 
66c328eeadSBarry Smith   /* Count currently needed memory usage and future memory requirements
67c328eeadSBarry Smith     (max, predicted)*/
689ccc4a9bSBarry Smith   for (i = i_row; i < nrows; i++) {
699ccc4a9bSBarry Smith     if (!result->row_nnz[i]) {
70b2f1ab58SBarry Smith       max_need_extra += max_row_nnz[i];
719ccc4a9bSBarry Smith     } else {
72b2f1ab58SBarry Smith       need_already += max_row_nnz[i];
73b2f1ab58SBarry Smith     }
74b2f1ab58SBarry Smith   }
75b2f1ab58SBarry Smith 
76c328eeadSBarry Smith   /* Make maximal and realistic memory requirement estimates */
77b2f1ab58SBarry Smith   n_alloc_max = n_alloc_ok + need_already + max_need_extra;
789ccc4a9bSBarry Smith   n_alloc_est = n_alloc_ok + need_already + (int)(((PetscReal)max_need_extra) * ((PetscReal)n_alloc_ok) / ((PetscReal)n_alloc_ok_max));
79b2f1ab58SBarry Smith 
80c328eeadSBarry Smith   /* Choose array sizes */
812205254eSKarl Rupp   if (n_alloc_max == n_alloc_est) n_alloc = n_alloc_max;
822205254eSKarl Rupp   else if (n_alloc_now >= n_alloc_est) n_alloc = n_alloc_now;
832205254eSKarl Rupp   else if (n_alloc_max < n_alloc_est * (1 + xtra_perc / 100.0)) n_alloc = n_alloc_max;
842205254eSKarl Rupp   else n_alloc = (int)(n_alloc_est * (1 + xtra_perc / 100.0));
85b2f1ab58SBarry Smith 
86c328eeadSBarry Smith   /* If new estimate is less than what we already have,
87c328eeadSBarry Smith     don't reallocate, just garbage-collect */
88ad540459SPierre Jolivet   if (n_alloc_max != n_alloc_est && n_alloc < result->n_alloc_icol) n_alloc = result->n_alloc_icol;
89b2f1ab58SBarry Smith 
90c328eeadSBarry Smith   /* Motivate dimension choice */
919566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "   Allocating %" PetscInt_FMT " nonzeros: ", n_alloc));
922205254eSKarl Rupp   if (n_alloc_max == n_alloc_est) {
939566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "this is the correct size\n"));
942205254eSKarl Rupp   } else if (n_alloc_now >= n_alloc_est) {
959566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "the current size, which seems enough\n"));
962205254eSKarl Rupp   } else if (n_alloc_max < n_alloc_est * (1 + xtra_perc / 100.0)) {
979566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "the maximum estimate\n"));
982205254eSKarl Rupp   } else {
999566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "%6.2f %% more than the estimate\n", (double)xtra_perc));
1002205254eSKarl Rupp   }
101b2f1ab58SBarry Smith 
102c328eeadSBarry Smith   /**********************************************************
103c328eeadSBarry Smith   2. Rescue arrays which would be lost
104c328eeadSBarry Smith   Count how many rows and nonzeros will have to be rescued
105c328eeadSBarry Smith   when moving all arrays in place */
106c599c493SJunchao Zhang   n_rescue = 0;
1072205254eSKarl Rupp   if (*n_row_alloc_ok == 0) *n_alloc_used = 0;
1089ccc4a9bSBarry Smith   else {
109b2f1ab58SBarry Smith     i = *n_row_alloc_ok - 1;
1102205254eSKarl Rupp 
1119ccc4a9bSBarry Smith     *n_alloc_used = (result->icols[i] - result->alloc_icol) + result->row_nnz[i];
112b2f1ab58SBarry Smith   }
113b2f1ab58SBarry Smith 
1149ccc4a9bSBarry Smith   for (i = *n_row_alloc_ok; i < nrows; i++) {
115b2f1ab58SBarry Smith     i_here = result->icols[i] - result->alloc_icol;
116b2f1ab58SBarry Smith     i_last = i_here + result->row_nnz[i];
1179ccc4a9bSBarry Smith     if (result->row_nnz[i] > 0) {
118ad540459SPierre Jolivet       if (*n_alloc_used > i_here || i_last > n_alloc) n_rescue += result->row_nnz[i];
119b2f1ab58SBarry Smith 
1202205254eSKarl Rupp       if (i < i_row) *n_alloc_used += result->row_nnz[i];
1212205254eSKarl Rupp       else *n_alloc_used += max_row_nnz[i];
122b2f1ab58SBarry Smith     }
123b2f1ab58SBarry Smith   }
124b2f1ab58SBarry Smith 
125c328eeadSBarry Smith   /* Allocate rescue arrays */
1269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_rescue, &icol_rescue));
1279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_rescue, &val_rescue));
128b2f1ab58SBarry Smith 
129c328eeadSBarry Smith   /* Rescue the arrays which need rescuing */
130c599c493SJunchao Zhang   n_rescue = 0;
1312205254eSKarl Rupp   if (*n_row_alloc_ok == 0) *n_alloc_used = 0;
1329ccc4a9bSBarry Smith   else {
133b2f1ab58SBarry Smith     i             = *n_row_alloc_ok - 1;
1349ccc4a9bSBarry Smith     *n_alloc_used = (result->icols[i] - result->alloc_icol) + result->row_nnz[i];
135b2f1ab58SBarry Smith   }
136b2f1ab58SBarry Smith 
1379ccc4a9bSBarry Smith   for (i = *n_row_alloc_ok; i < nrows; i++) {
138b2f1ab58SBarry Smith     i_here = result->icols[i] - result->alloc_icol;
139b2f1ab58SBarry Smith     i_last = i_here + result->row_nnz[i];
1409ccc4a9bSBarry Smith     if (result->row_nnz[i] > 0) {
1414e1ff37aSBarry Smith       if (*n_alloc_used > i_here || i_last > n_alloc) {
1429566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&icol_rescue[n_rescue], result->icols[i], result->row_nnz[i]));
1439566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&val_rescue[n_rescue], result->values[i], result->row_nnz[i]));
144b2f1ab58SBarry Smith         n_rescue += result->row_nnz[i];
145b2f1ab58SBarry Smith       }
146b2f1ab58SBarry Smith 
1472205254eSKarl Rupp       if (i < i_row) *n_alloc_used += result->row_nnz[i];
1482205254eSKarl Rupp       else *n_alloc_used += max_row_nnz[i];
149b2f1ab58SBarry Smith     }
150b2f1ab58SBarry Smith   }
151b2f1ab58SBarry Smith 
152c328eeadSBarry Smith   /**********************************************************
153c328eeadSBarry Smith   3. Reallocate and correct administration */
154b2f1ab58SBarry Smith 
1559ccc4a9bSBarry Smith   if (n_alloc != result->n_alloc_icol) {
156cc6fc633SBarry Smith     n_copy = PetscMin(n_alloc, result->n_alloc_icol);
157b2f1ab58SBarry Smith 
158c328eeadSBarry Smith     /* PETSC knows no REALLOC, so we'll REALLOC ourselves.
159b2f1ab58SBarry Smith 
160c328eeadSBarry Smith         Allocate new icol-data, copy old contents */
1619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n_alloc, &result->alloc_icol));
1629566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(result->alloc_icol, alloc_icol_old, n_copy));
163b2f1ab58SBarry Smith 
164c328eeadSBarry Smith     /* Update administration, Reset pointers to new arrays  */
165b2f1ab58SBarry Smith     result->n_alloc_icol = n_alloc;
1669ccc4a9bSBarry Smith     for (i = 0; i < nrows; i++) {
1679ccc4a9bSBarry Smith       result->icols[i]  = result->alloc_icol + (result->icols[i] - alloc_icol_old);
1689ccc4a9bSBarry Smith       result->values[i] = result->alloc_val + (result->icols[i] - result->alloc_icol);
169b2f1ab58SBarry Smith     }
170b2f1ab58SBarry Smith 
171c328eeadSBarry Smith     /* Delete old array */
1729566063dSJacob Faibussowitsch     PetscCall(PetscFree(alloc_icol_old));
173b2f1ab58SBarry Smith 
174c328eeadSBarry Smith     /* Allocate new value-data, copy old contents */
1759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n_alloc, &result->alloc_val));
1769566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(result->alloc_val, alloc_val_old, n_copy));
177b2f1ab58SBarry Smith 
178c328eeadSBarry Smith     /* Update administration, Reset pointers to new arrays  */
179b2f1ab58SBarry Smith     result->n_alloc_val = n_alloc;
180ad540459SPierre Jolivet     for (i = 0; i < nrows; i++) result->values[i] = result->alloc_val + (result->icols[i] - result->alloc_icol);
181b2f1ab58SBarry Smith 
182c328eeadSBarry Smith     /* Delete old array */
1839566063dSJacob Faibussowitsch     PetscCall(PetscFree(alloc_val_old));
184b2f1ab58SBarry Smith   }
185b2f1ab58SBarry Smith 
186c328eeadSBarry Smith   /*********************************************************
187c328eeadSBarry Smith   4. Copy all the arrays to their proper places */
188c599c493SJunchao Zhang   n_rescue = 0;
1892205254eSKarl Rupp   if (*n_row_alloc_ok == 0) *n_alloc_used = 0;
1909ccc4a9bSBarry Smith   else {
191b2f1ab58SBarry Smith     i = *n_row_alloc_ok - 1;
1922205254eSKarl Rupp 
1939ccc4a9bSBarry Smith     *n_alloc_used = (result->icols[i] - result->alloc_icol) + result->row_nnz[i];
194b2f1ab58SBarry Smith   }
195b2f1ab58SBarry Smith 
1969ccc4a9bSBarry Smith   for (i = *n_row_alloc_ok; i < nrows; i++) {
197b2f1ab58SBarry Smith     i_here = result->icols[i] - result->alloc_icol;
198b2f1ab58SBarry Smith     i_last = i_here + result->row_nnz[i];
199b2f1ab58SBarry Smith 
200b2f1ab58SBarry Smith     result->icols[i]  = result->alloc_icol + *n_alloc_used;
201b2f1ab58SBarry Smith     result->values[i] = result->alloc_val + *n_alloc_used;
202b2f1ab58SBarry Smith 
2039ccc4a9bSBarry Smith     if (result->row_nnz[i] > 0) {
2049ccc4a9bSBarry Smith       if (*n_alloc_used > i_here || i_last > n_alloc) {
2059566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(result->icols[i], &icol_rescue[n_rescue], result->row_nnz[i]));
2069566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(result->values[i], &val_rescue[n_rescue], result->row_nnz[i]));
2072205254eSKarl Rupp 
208b2f1ab58SBarry Smith         n_rescue += result->row_nnz[i];
2099ccc4a9bSBarry Smith       } else {
2109ccc4a9bSBarry Smith         for (j = 0; j < result->row_nnz[i]; j++) {
211b2f1ab58SBarry Smith           result->icols[i][j]  = result->alloc_icol[i_here + j];
212b2f1ab58SBarry Smith           result->values[i][j] = result->alloc_val[i_here + j];
213b2f1ab58SBarry Smith         }
214b2f1ab58SBarry Smith       }
2152205254eSKarl Rupp       if (i < i_row) *n_alloc_used += result->row_nnz[i];
2162205254eSKarl Rupp       else *n_alloc_used += max_row_nnz[i];
217b2f1ab58SBarry Smith     }
218b2f1ab58SBarry Smith   }
219b2f1ab58SBarry Smith 
220c328eeadSBarry Smith   /* Delete the rescue arrays */
2219566063dSJacob Faibussowitsch   PetscCall(PetscFree(icol_rescue));
2229566063dSJacob Faibussowitsch   PetscCall(PetscFree(val_rescue));
2232205254eSKarl Rupp 
224b2f1ab58SBarry Smith   *n_row_alloc_ok = i_row;
225b2f1ab58SBarry Smith   PetscFunctionReturn(0);
226b2f1ab58SBarry Smith }
227b2f1ab58SBarry Smith 
228b2f1ab58SBarry Smith /*
229b2f1ab58SBarry Smith   spbas_incomplete_cholesky:
230b2f1ab58SBarry Smith      incomplete Cholesky decomposition of a square, symmetric,
231b2f1ab58SBarry Smith      positive definite matrix.
232b2f1ab58SBarry Smith 
233b2f1ab58SBarry Smith      In case negative diagonals are encountered, function returns
234b2f1ab58SBarry Smith      NEGATIVE_DIAGONAL. When this happens, call this function again
235b2f1ab58SBarry Smith      with a larger epsdiag_in, a less sparse pattern, and/or a smaller
236b2f1ab58SBarry Smith      droptol
237b2f1ab58SBarry Smith */
2389371c9d4SSatish Balay 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) {
239b2f1ab58SBarry Smith   PetscInt        jL;
240b2f1ab58SBarry Smith   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *)A->data;
241b2f1ab58SBarry Smith   PetscInt       *ai = a->i, *aj = a->j;
242b2f1ab58SBarry Smith   MatScalar      *aa = a->a;
243b2f1ab58SBarry Smith   PetscInt        nrows, ncols;
244b2f1ab58SBarry Smith   PetscInt       *max_row_nnz;
245b2f1ab58SBarry Smith   spbas_matrix    retval;
246b2f1ab58SBarry Smith   PetscScalar    *diag;
247b2f1ab58SBarry Smith   PetscScalar    *val;
248b2f1ab58SBarry Smith   PetscScalar    *lvec;
249b2f1ab58SBarry Smith   PetscScalar     epsdiag;
250b2f1ab58SBarry Smith   PetscInt        i, j, k;
251ace3abfcSBarry Smith   const PetscBool do_values = PETSC_TRUE;
2529ccc4a9bSBarry Smith   PetscInt       *r1_icol;
2539ccc4a9bSBarry Smith   PetscScalar    *r1_val;
2549ccc4a9bSBarry Smith   PetscInt       *r_icol;
2559ccc4a9bSBarry Smith   PetscInt        r_nnz;
2569ccc4a9bSBarry Smith   PetscScalar    *r_val;
2579ccc4a9bSBarry Smith   PetscInt       *A_icol;
2589ccc4a9bSBarry Smith   PetscInt        A_nnz;
2599ccc4a9bSBarry Smith   PetscScalar    *A_val;
2609ccc4a9bSBarry Smith   PetscInt       *p_icol;
2619ccc4a9bSBarry Smith   PetscInt        p_nnz;
2629ccc4a9bSBarry Smith   PetscInt        n_row_alloc_ok = 0; /* number of rows which have been stored   correctly in the matrix */
2639ccc4a9bSBarry Smith   PetscInt        n_alloc_used   = 0; /* part of result->icols and result->values   which is currently being used */
264b2f1ab58SBarry Smith 
2659ccc4a9bSBarry Smith   PetscFunctionBegin;
2664e1ff37aSBarry Smith   /* Convert the Manteuffel shift from 'fraction of average diagonal' to   dimensioned value */
2679566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &nrows, &ncols));
2689566063dSJacob Faibussowitsch   PetscCall(MatGetTrace(A, &epsdiag));
2692205254eSKarl Rupp 
270b2f1ab58SBarry Smith   epsdiag *= epsdiag_in / nrows;
2712205254eSKarl Rupp 
2729566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "   Dimensioned Manteuffel shift %g Drop tolerance %g\n", (double)PetscRealPart(epsdiag), (double)droptol));
273b2f1ab58SBarry Smith 
2742205254eSKarl Rupp   if (droptol < 1e-10) droptol = 1e-10;
275b2f1ab58SBarry Smith 
276b2f1ab58SBarry Smith   retval.nrows        = nrows;
277b2f1ab58SBarry Smith   retval.ncols        = nrows;
278b2f1ab58SBarry Smith   retval.nnz          = pattern.nnz / 10;
279b2f1ab58SBarry Smith   retval.col_idx_type = SPBAS_COLUMN_NUMBERS;
280b2f1ab58SBarry Smith   retval.block_data   = PETSC_TRUE;
281b2f1ab58SBarry Smith 
2829566063dSJacob Faibussowitsch   PetscCall(spbas_allocate_pattern(&retval, do_values));
2839566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(retval.row_nnz, nrows));
2849566063dSJacob Faibussowitsch   PetscCall(spbas_allocate_data(&retval));
285b2f1ab58SBarry Smith   retval.nnz = 0;
286b2f1ab58SBarry Smith 
2879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrows, &diag));
2889566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nrows, &val));
2899566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nrows, &lvec));
2909566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nrows, &max_row_nnz));
291b2f1ab58SBarry Smith 
292c328eeadSBarry Smith   /* Count the nonzeros on transpose of pattern */
2934e1ff37aSBarry Smith   for (i = 0; i < nrows; i++) {
294b2f1ab58SBarry Smith     p_nnz  = pattern.row_nnz[i];
295b2f1ab58SBarry Smith     p_icol = pattern.icols[i];
296ad540459SPierre Jolivet     for (j = 0; j < p_nnz; j++) max_row_nnz[i + p_icol[j]]++;
297b2f1ab58SBarry Smith   }
298b2f1ab58SBarry Smith 
299c328eeadSBarry Smith   /* Calculate rows of L */
3004e1ff37aSBarry Smith   for (i = 0; i < nrows; i++) {
301b2f1ab58SBarry Smith     p_nnz  = pattern.row_nnz[i];
302b2f1ab58SBarry Smith     p_icol = pattern.icols[i];
303b2f1ab58SBarry Smith 
304b2f1ab58SBarry Smith     r_nnz  = retval.row_nnz[i];
305b2f1ab58SBarry Smith     r_icol = retval.icols[i];
306b2f1ab58SBarry Smith     r_val  = retval.values[i];
307b2f1ab58SBarry Smith     A_nnz  = ai[rip[i] + 1] - ai[rip[i]];
308b2f1ab58SBarry Smith     A_icol = &aj[ai[rip[i]]];
309b2f1ab58SBarry Smith     A_val  = &aa[ai[rip[i]]];
310b2f1ab58SBarry Smith 
311c328eeadSBarry Smith     /* Calculate  val += A(i,i:n)'; */
3124e1ff37aSBarry Smith     for (j = 0; j < A_nnz; j++) {
313b2f1ab58SBarry Smith       k = riip[A_icol[j]];
3142205254eSKarl Rupp       if (k >= i) val[k] = A_val[j];
315b2f1ab58SBarry Smith     }
316b2f1ab58SBarry Smith 
317c328eeadSBarry Smith     /*  Add regularization */
318b2f1ab58SBarry Smith     val[i] += epsdiag;
319b2f1ab58SBarry Smith 
320c328eeadSBarry Smith     /* Calculate lvec   = diag(D(0:i-1)) * L(0:i-1,i);
321c328eeadSBarry Smith         val(i) = A(i,i) - L(0:i-1,i)' * lvec */
3224e1ff37aSBarry Smith     for (j = 0; j < r_nnz; j++) {
323b2f1ab58SBarry Smith       k       = r_icol[j];
324b2f1ab58SBarry Smith       lvec[k] = diag[k] * r_val[j];
325b2f1ab58SBarry Smith       val[i] -= r_val[j] * lvec[k];
326b2f1ab58SBarry Smith     }
327b2f1ab58SBarry Smith 
328c328eeadSBarry Smith     /* Calculate the new diagonal */
329b2f1ab58SBarry Smith     diag[i] = val[i];
3304e1ff37aSBarry Smith     if (PetscRealPart(diag[i]) < droptol) {
3319566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Error in spbas_incomplete_cholesky:\n"));
3329566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Negative diagonal in row %" PetscInt_FMT "\n", i + 1));
333b2f1ab58SBarry Smith 
334c328eeadSBarry Smith       /* Delete the whole matrix at once. */
3359566063dSJacob Faibussowitsch       PetscCall(spbas_delete(retval));
336cc442ecdSBarry Smith       *success = PETSC_FALSE;
337cc442ecdSBarry Smith       PetscFunctionReturn(0);
338b2f1ab58SBarry Smith     }
339b2f1ab58SBarry Smith 
340c328eeadSBarry Smith     /* If necessary, allocate arrays */
3414e1ff37aSBarry Smith     if (r_nnz == 0) {
342cc442ecdSBarry Smith       PetscBool success = spbas_cholesky_row_alloc(retval, i, 1, &n_alloc_used);
34328b400f6SJacob Faibussowitsch       PetscCheck(success, PETSC_COMM_SELF, PETSC_ERR_MEM, "spbas_cholesky_row_alloc() failed");
344b2f1ab58SBarry Smith       r_icol = retval.icols[i];
345b2f1ab58SBarry Smith       r_val  = retval.values[i];
346b2f1ab58SBarry Smith     }
347b2f1ab58SBarry Smith 
348c328eeadSBarry Smith     /* Now, fill in */
349b2f1ab58SBarry Smith     r_icol[r_nnz] = i;
350b2f1ab58SBarry Smith     r_val[r_nnz]  = 1.0;
351b2f1ab58SBarry Smith     r_nnz++;
352b2f1ab58SBarry Smith     retval.row_nnz[i]++;
353b2f1ab58SBarry Smith 
354b2f1ab58SBarry Smith     retval.nnz += r_nnz;
355b2f1ab58SBarry Smith 
356c328eeadSBarry Smith     /* Calculate
357c328eeadSBarry Smith         val(i+1:n) = (A(i,i+1:n)- L(0:i-1,i+1:n)' * lvec)/diag(i) */
3584e1ff37aSBarry Smith     for (j = 1; j < p_nnz; j++) {
359b2f1ab58SBarry Smith       k       = i + p_icol[j];
360b2f1ab58SBarry Smith       r1_icol = retval.icols[k];
361b2f1ab58SBarry Smith       r1_val  = retval.values[k];
362ad540459SPierre Jolivet       for (jL = 0; jL < retval.row_nnz[k]; jL++) val[k] -= r1_val[jL] * lvec[r1_icol[jL]];
363b2f1ab58SBarry Smith       val[k] /= diag[i];
364b2f1ab58SBarry Smith 
3654e1ff37aSBarry Smith       if (PetscAbsScalar(val[k]) > droptol || PetscAbsScalar(val[k]) < -droptol) {
366c328eeadSBarry Smith         /* If necessary, allocate arrays */
367cc442ecdSBarry Smith         if (!retval.row_nnz[k]) {
368cc442ecdSBarry Smith           PetscBool flag, success = spbas_cholesky_row_alloc(retval, k, max_row_nnz[k], &n_alloc_used);
369cc442ecdSBarry Smith           if (!success) {
3709566063dSJacob Faibussowitsch             PetscCall(spbas_cholesky_garbage_collect(&retval, i, &n_row_alloc_ok, &n_alloc_used, max_row_nnz));
371cc442ecdSBarry Smith             flag = spbas_cholesky_row_alloc(retval, k, max_row_nnz[k], &n_alloc_used);
37228b400f6SJacob Faibussowitsch             PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_MEM, "Allocation in spbas_cholesky_row_alloc() failed");
373b2f1ab58SBarry Smith             r_icol = retval.icols[i];
374cc442ecdSBarry Smith           }
375b2f1ab58SBarry Smith         }
376b2f1ab58SBarry Smith 
377b2f1ab58SBarry Smith         retval.icols[k][retval.row_nnz[k]]  = i;
378b2f1ab58SBarry Smith         retval.values[k][retval.row_nnz[k]] = val[k];
379b2f1ab58SBarry Smith         retval.row_nnz[k]++;
380b2f1ab58SBarry Smith       }
381b2f1ab58SBarry Smith       val[k] = 0;
382b2f1ab58SBarry Smith     }
383b2f1ab58SBarry Smith 
38471d55d18SBarry Smith     /* Erase the values used in the work arrays */
3852205254eSKarl Rupp     for (j = 0; j < r_nnz; j++) lvec[r_icol[j]] = 0;
386b2f1ab58SBarry Smith   }
387b2f1ab58SBarry Smith 
3889566063dSJacob Faibussowitsch   PetscCall(PetscFree(lvec));
3899566063dSJacob Faibussowitsch   PetscCall(PetscFree(val));
390b2f1ab58SBarry Smith 
3919566063dSJacob Faibussowitsch   PetscCall(spbas_cholesky_garbage_collect(&retval, nrows, &n_row_alloc_ok, &n_alloc_used, max_row_nnz));
3929566063dSJacob Faibussowitsch   PetscCall(PetscFree(max_row_nnz));
393b2f1ab58SBarry Smith 
394c328eeadSBarry Smith   /* Place the inverse of the diagonals in the matrix */
3959ccc4a9bSBarry Smith   for (i = 0; i < nrows; i++) {
396b2f1ab58SBarry Smith     r_nnz = retval.row_nnz[i];
3972205254eSKarl Rupp 
398b2f1ab58SBarry Smith     retval.values[i][r_nnz - 1] = 1.0 / diag[i];
399ad540459SPierre Jolivet     for (j = 0; j < r_nnz - 1; j++) retval.values[i][j] *= -1;
400b2f1ab58SBarry Smith   }
4019566063dSJacob Faibussowitsch   PetscCall(PetscFree(diag));
402b2f1ab58SBarry Smith   *matrix_L = retval;
403cc442ecdSBarry Smith   *success  = PETSC_TRUE;
404b2f1ab58SBarry Smith   PetscFunctionReturn(0);
405b2f1ab58SBarry Smith }
406