xref: /petsc/src/mat/impls/aij/seq/bas/spbas_cholesky.h (revision c599c4937b1a20dc44a55edd8eb2a391549cfdf4)
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        max_need_extra= 0;
50b2f1ab58SBarry Smith   PetscInt        n_alloc_max, n_alloc_est, n_alloc;
51b2f1ab58SBarry Smith   PetscInt        n_alloc_now     = result->n_alloc_icol;
52b2f1ab58SBarry Smith   PetscInt        *alloc_icol_old = result->alloc_icol;
53b2f1ab58SBarry Smith   PetscScalar     *alloc_val_old  = result->alloc_val;
54b2f1ab58SBarry Smith   PetscInt        *icol_rescue;
55b2f1ab58SBarry Smith   PetscScalar     *val_rescue;
56b2f1ab58SBarry Smith   PetscInt        n_rescue;
57b2f1ab58SBarry Smith   PetscInt        i_here, i_last, n_copy;
58d36aa34eSBarry Smith   const PetscReal xtra_perc = 20;
59b2f1ab58SBarry Smith 
60b2f1ab58SBarry Smith   PetscFunctionBegin;
61c328eeadSBarry Smith   /*********************************************************
62c328eeadSBarry Smith   1. Choose appropriate array size
63c328eeadSBarry Smith   Count number of rows and memory usage which is already final */
649ccc4a9bSBarry Smith   for (i=0; i<i_row; i++)  {
65b2f1ab58SBarry Smith     n_alloc_ok     += result->row_nnz[i];
66b2f1ab58SBarry Smith     n_alloc_ok_max += max_row_nnz[i];
67b2f1ab58SBarry Smith   }
68b2f1ab58SBarry Smith 
69c328eeadSBarry Smith   /* Count currently needed memory usage and future memory requirements
70c328eeadSBarry Smith     (max, predicted)*/
719ccc4a9bSBarry Smith   for (i=i_row; i<nrows; i++) {
729ccc4a9bSBarry Smith     if (!result->row_nnz[i]) {
73b2f1ab58SBarry Smith       max_need_extra += max_row_nnz[i];
749ccc4a9bSBarry Smith     } else {
75b2f1ab58SBarry Smith       need_already += max_row_nnz[i];
76b2f1ab58SBarry Smith     }
77b2f1ab58SBarry Smith   }
78b2f1ab58SBarry Smith 
79c328eeadSBarry Smith   /* Make maximal and realistic memory requirement estimates */
80b2f1ab58SBarry Smith   n_alloc_max = n_alloc_ok + need_already + max_need_extra;
819ccc4a9bSBarry Smith   n_alloc_est = n_alloc_ok + need_already + (int) (((PetscReal) max_need_extra) *  ((PetscReal) n_alloc_ok) /((PetscReal) n_alloc_ok_max));
82b2f1ab58SBarry Smith 
83c328eeadSBarry Smith   /* Choose array sizes */
842205254eSKarl Rupp   if (n_alloc_max == n_alloc_est) n_alloc = n_alloc_max;
852205254eSKarl Rupp   else if (n_alloc_now >= n_alloc_est) n_alloc = n_alloc_now;
862205254eSKarl Rupp   else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0)) n_alloc = n_alloc_max;
872205254eSKarl Rupp   else n_alloc = (int) (n_alloc_est * (1+xtra_perc/100.0));
88b2f1ab58SBarry Smith 
89c328eeadSBarry Smith   /* If new estimate is less than what we already have,
90c328eeadSBarry Smith     don't reallocate, just garbage-collect */
919ccc4a9bSBarry Smith   if (n_alloc_max != n_alloc_est && n_alloc < result->n_alloc_icol) {
92b2f1ab58SBarry Smith     n_alloc = result->n_alloc_icol;
93b2f1ab58SBarry Smith   }
94b2f1ab58SBarry Smith 
95c328eeadSBarry Smith   /* Motivate dimension choice */
969566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL,"   Allocating %" PetscInt_FMT " nonzeros: ",n_alloc));
972205254eSKarl Rupp   if (n_alloc_max == n_alloc_est) {
989566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL,"this is the correct size\n"));
992205254eSKarl Rupp   } else if (n_alloc_now >= n_alloc_est) {
1009566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL,"the current size, which seems enough\n"));
1012205254eSKarl Rupp   } else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0)) {
1029566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL,"the maximum estimate\n"));
1032205254eSKarl Rupp   } else {
1049566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL,"%6.2f %% more than the estimate\n",(double)xtra_perc));
1052205254eSKarl Rupp   }
106b2f1ab58SBarry Smith 
107c328eeadSBarry Smith   /**********************************************************
108c328eeadSBarry Smith   2. Rescue arrays which would be lost
109c328eeadSBarry Smith   Count how many rows and nonzeros will have to be rescued
110c328eeadSBarry Smith   when moving all arrays in place */
111*c599c493SJunchao Zhang   n_rescue = 0;
1122205254eSKarl Rupp   if (*n_row_alloc_ok==0) *n_alloc_used = 0;
1139ccc4a9bSBarry Smith   else {
114b2f1ab58SBarry Smith     i = *n_row_alloc_ok - 1;
1152205254eSKarl Rupp 
1169ccc4a9bSBarry Smith     *n_alloc_used = (result->icols[i]-result->alloc_icol) +  result->row_nnz[i];
117b2f1ab58SBarry Smith   }
118b2f1ab58SBarry Smith 
1199ccc4a9bSBarry Smith   for (i=*n_row_alloc_ok; i<nrows; i++) {
120b2f1ab58SBarry Smith     i_here = result->icols[i]-result->alloc_icol;
121b2f1ab58SBarry Smith     i_last = i_here + result->row_nnz[i];
1229ccc4a9bSBarry Smith     if (result->row_nnz[i]>0) {
1239ccc4a9bSBarry Smith       if (*n_alloc_used > i_here || i_last > n_alloc) {
124b2f1ab58SBarry Smith         n_rescue += result->row_nnz[i];
125b2f1ab58SBarry Smith       }
126b2f1ab58SBarry Smith 
1272205254eSKarl Rupp       if (i<i_row) *n_alloc_used += result->row_nnz[i];
1282205254eSKarl Rupp       else         *n_alloc_used += max_row_nnz[i];
129b2f1ab58SBarry Smith     }
130b2f1ab58SBarry Smith   }
131b2f1ab58SBarry Smith 
132c328eeadSBarry Smith   /* Allocate rescue arrays */
1339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_rescue, &icol_rescue));
1349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_rescue, &val_rescue));
135b2f1ab58SBarry Smith 
136c328eeadSBarry Smith   /* Rescue the arrays which need rescuing */
137*c599c493SJunchao Zhang   n_rescue = 0;
1382205254eSKarl Rupp   if (*n_row_alloc_ok==0) *n_alloc_used = 0;
1399ccc4a9bSBarry Smith   else {
140b2f1ab58SBarry Smith     i             = *n_row_alloc_ok - 1;
1419ccc4a9bSBarry Smith     *n_alloc_used = (result->icols[i]-result->alloc_icol) +  result->row_nnz[i];
142b2f1ab58SBarry Smith   }
143b2f1ab58SBarry Smith 
1449ccc4a9bSBarry Smith   for (i=*n_row_alloc_ok; i<nrows; i++) {
145b2f1ab58SBarry Smith     i_here = result->icols[i]-result->alloc_icol;
146b2f1ab58SBarry Smith     i_last = i_here + result->row_nnz[i];
1479ccc4a9bSBarry Smith     if (result->row_nnz[i]>0) {
1484e1ff37aSBarry Smith       if (*n_alloc_used > i_here || i_last > n_alloc) {
1499566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&icol_rescue[n_rescue], result->icols[i], result->row_nnz[i]));
1509566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&val_rescue[n_rescue], result->values[i], result->row_nnz[i]));
151b2f1ab58SBarry Smith         n_rescue += result->row_nnz[i];
152b2f1ab58SBarry Smith       }
153b2f1ab58SBarry Smith 
1542205254eSKarl Rupp       if (i<i_row) *n_alloc_used += result->row_nnz[i];
1552205254eSKarl Rupp       else         *n_alloc_used += max_row_nnz[i];
156b2f1ab58SBarry Smith     }
157b2f1ab58SBarry Smith   }
158b2f1ab58SBarry Smith 
159c328eeadSBarry Smith   /**********************************************************
160c328eeadSBarry Smith   3. Reallocate and correct administration */
161b2f1ab58SBarry Smith 
1629ccc4a9bSBarry Smith   if (n_alloc != result->n_alloc_icol) {
163cc6fc633SBarry Smith     n_copy = PetscMin(n_alloc,result->n_alloc_icol);
164b2f1ab58SBarry Smith 
165c328eeadSBarry Smith     /* PETSC knows no REALLOC, so we'll REALLOC ourselves.
166b2f1ab58SBarry Smith 
167c328eeadSBarry Smith         Allocate new icol-data, copy old contents */
1689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n_alloc, &result->alloc_icol));
1699566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(result->alloc_icol, alloc_icol_old, n_copy));
170b2f1ab58SBarry Smith 
171c328eeadSBarry Smith     /* Update administration, Reset pointers to new arrays  */
172b2f1ab58SBarry Smith     result->n_alloc_icol = n_alloc;
1739ccc4a9bSBarry Smith     for (i=0; i<nrows; i++) {
1749ccc4a9bSBarry Smith       result->icols[i]  =  result->alloc_icol + (result->icols[i]  - alloc_icol_old);
1759ccc4a9bSBarry Smith       result->values[i] =  result->alloc_val  + (result->icols[i]  - result->alloc_icol);
176b2f1ab58SBarry Smith     }
177b2f1ab58SBarry Smith 
178c328eeadSBarry Smith     /* Delete old array */
1799566063dSJacob Faibussowitsch     PetscCall(PetscFree(alloc_icol_old));
180b2f1ab58SBarry Smith 
181c328eeadSBarry Smith     /* Allocate new value-data, copy old contents */
1829566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n_alloc, &result->alloc_val));
1839566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(result->alloc_val, alloc_val_old, n_copy));
184b2f1ab58SBarry Smith 
185c328eeadSBarry Smith     /* Update administration, Reset pointers to new arrays  */
186b2f1ab58SBarry Smith     result->n_alloc_val = n_alloc;
1879ccc4a9bSBarry Smith     for (i=0; i<nrows; i++) {
1889ccc4a9bSBarry Smith       result->values[i] =  result->alloc_val + (result->icols[i]  - result->alloc_icol);
189b2f1ab58SBarry Smith     }
190b2f1ab58SBarry Smith 
191c328eeadSBarry Smith     /* Delete old array */
1929566063dSJacob Faibussowitsch     PetscCall(PetscFree(alloc_val_old));
193b2f1ab58SBarry Smith   }
194b2f1ab58SBarry Smith 
195c328eeadSBarry Smith   /*********************************************************
196c328eeadSBarry Smith   4. Copy all the arrays to their proper places */
197*c599c493SJunchao Zhang   n_rescue = 0;
1982205254eSKarl Rupp   if (*n_row_alloc_ok==0) *n_alloc_used = 0;
1999ccc4a9bSBarry Smith   else {
200b2f1ab58SBarry Smith     i = *n_row_alloc_ok - 1;
2012205254eSKarl Rupp 
2029ccc4a9bSBarry Smith     *n_alloc_used = (result->icols[i]-result->alloc_icol) +  result->row_nnz[i];
203b2f1ab58SBarry Smith   }
204b2f1ab58SBarry Smith 
2059ccc4a9bSBarry Smith   for (i=*n_row_alloc_ok; i<nrows; i++) {
206b2f1ab58SBarry Smith     i_here = result->icols[i]-result->alloc_icol;
207b2f1ab58SBarry Smith     i_last = i_here + result->row_nnz[i];
208b2f1ab58SBarry Smith 
209b2f1ab58SBarry Smith     result->icols[i] = result->alloc_icol + *n_alloc_used;
210b2f1ab58SBarry Smith     result->values[i]= result->alloc_val  + *n_alloc_used;
211b2f1ab58SBarry Smith 
2129ccc4a9bSBarry Smith     if (result->row_nnz[i]>0) {
2139ccc4a9bSBarry Smith       if (*n_alloc_used > i_here || i_last > n_alloc) {
2149566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(result->icols[i], &icol_rescue[n_rescue], result->row_nnz[i]));
2159566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(result->values[i],&val_rescue[n_rescue],result->row_nnz[i]));
2162205254eSKarl Rupp 
217b2f1ab58SBarry Smith         n_rescue += result->row_nnz[i];
2189ccc4a9bSBarry Smith       } else {
2199ccc4a9bSBarry Smith         for (j=0; j<result->row_nnz[i]; j++) {
220b2f1ab58SBarry Smith           result->icols[i][j]  = result->alloc_icol[i_here+j];
221b2f1ab58SBarry Smith           result->values[i][j] = result->alloc_val[i_here+j];
222b2f1ab58SBarry Smith         }
223b2f1ab58SBarry Smith       }
2242205254eSKarl Rupp       if (i<i_row) *n_alloc_used += result->row_nnz[i];
2252205254eSKarl Rupp       else         *n_alloc_used += max_row_nnz[i];
226b2f1ab58SBarry Smith     }
227b2f1ab58SBarry Smith   }
228b2f1ab58SBarry Smith 
229c328eeadSBarry Smith   /* Delete the rescue arrays */
2309566063dSJacob Faibussowitsch   PetscCall(PetscFree(icol_rescue));
2319566063dSJacob Faibussowitsch   PetscCall(PetscFree(val_rescue));
2322205254eSKarl Rupp 
233b2f1ab58SBarry Smith   *n_row_alloc_ok = i_row;
234b2f1ab58SBarry Smith   PetscFunctionReturn(0);
235b2f1ab58SBarry Smith }
236b2f1ab58SBarry Smith 
237b2f1ab58SBarry Smith /*
238b2f1ab58SBarry Smith   spbas_incomplete_cholesky:
239b2f1ab58SBarry Smith      incomplete Cholesky decomposition of a square, symmetric,
240b2f1ab58SBarry Smith      positive definite matrix.
241b2f1ab58SBarry Smith 
242b2f1ab58SBarry Smith      In case negative diagonals are encountered, function returns
243b2f1ab58SBarry Smith      NEGATIVE_DIAGONAL. When this happens, call this function again
244b2f1ab58SBarry Smith      with a larger epsdiag_in, a less sparse pattern, and/or a smaller
245b2f1ab58SBarry Smith      droptol
246b2f1ab58SBarry Smith */
247cc442ecdSBarry 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)
248b2f1ab58SBarry Smith {
249b2f1ab58SBarry Smith   PetscInt        jL;
250b2f1ab58SBarry Smith   Mat_SeqAIJ      *a =(Mat_SeqAIJ*)A->data;
251b2f1ab58SBarry Smith   PetscInt        *ai=a->i,*aj=a->j;
252b2f1ab58SBarry Smith   MatScalar       *aa=a->a;
253b2f1ab58SBarry Smith   PetscInt        nrows, ncols;
254b2f1ab58SBarry Smith   PetscInt        *max_row_nnz;
255b2f1ab58SBarry Smith   spbas_matrix    retval;
256b2f1ab58SBarry Smith   PetscScalar     *diag;
257b2f1ab58SBarry Smith   PetscScalar     *val;
258b2f1ab58SBarry Smith   PetscScalar     *lvec;
259b2f1ab58SBarry Smith   PetscScalar     epsdiag;
260b2f1ab58SBarry Smith   PetscInt        i,j,k;
261ace3abfcSBarry Smith   const PetscBool do_values = PETSC_TRUE;
2629ccc4a9bSBarry Smith   PetscInt        *r1_icol;
2639ccc4a9bSBarry Smith   PetscScalar     *r1_val;
2649ccc4a9bSBarry Smith   PetscInt        *r_icol;
2659ccc4a9bSBarry Smith   PetscInt        r_nnz;
2669ccc4a9bSBarry Smith   PetscScalar     *r_val;
2679ccc4a9bSBarry Smith   PetscInt        *A_icol;
2689ccc4a9bSBarry Smith   PetscInt        A_nnz;
2699ccc4a9bSBarry Smith   PetscScalar     *A_val;
2709ccc4a9bSBarry Smith   PetscInt        *p_icol;
2719ccc4a9bSBarry Smith   PetscInt        p_nnz;
2729ccc4a9bSBarry Smith   PetscInt        n_row_alloc_ok = 0;   /* number of rows which have been stored   correctly in the matrix */
2739ccc4a9bSBarry Smith   PetscInt        n_alloc_used   = 0;   /* part of result->icols and result->values   which is currently being used */
274b2f1ab58SBarry Smith 
2759ccc4a9bSBarry Smith   PetscFunctionBegin;
2764e1ff37aSBarry Smith   /* Convert the Manteuffel shift from 'fraction of average diagonal' to   dimensioned value */
2779566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &nrows, &ncols));
2789566063dSJacob Faibussowitsch   PetscCall(MatGetTrace(A, &epsdiag));
2792205254eSKarl Rupp 
280b2f1ab58SBarry Smith   epsdiag *= epsdiag_in / nrows;
2812205254eSKarl Rupp 
2829566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL,"   Dimensioned Manteuffel shift %g Drop tolerance %g\n", (double)PetscRealPart(epsdiag),(double)droptol));
283b2f1ab58SBarry Smith 
2842205254eSKarl Rupp   if (droptol<1e-10) droptol=1e-10;
285b2f1ab58SBarry Smith 
286b2f1ab58SBarry Smith   retval.nrows        = nrows;
287b2f1ab58SBarry Smith   retval.ncols        = nrows;
288b2f1ab58SBarry Smith   retval.nnz          = pattern.nnz/10;
289b2f1ab58SBarry Smith   retval.col_idx_type = SPBAS_COLUMN_NUMBERS;
290b2f1ab58SBarry Smith   retval.block_data   = PETSC_TRUE;
291b2f1ab58SBarry Smith 
2929566063dSJacob Faibussowitsch   PetscCall(spbas_allocate_pattern(&retval, do_values));
2939566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(retval.row_nnz, nrows));
2949566063dSJacob Faibussowitsch   PetscCall(spbas_allocate_data(&retval));
295b2f1ab58SBarry Smith   retval.nnz = 0;
296b2f1ab58SBarry Smith 
2979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrows, &diag));
2989566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nrows, &val));
2999566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nrows, &lvec));
3009566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nrows, &max_row_nnz));
301b2f1ab58SBarry Smith 
302c328eeadSBarry Smith   /* Count the nonzeros on transpose of pattern */
3034e1ff37aSBarry Smith   for (i = 0; i<nrows; i++)  {
304b2f1ab58SBarry Smith     p_nnz  = pattern.row_nnz[i];
305b2f1ab58SBarry Smith     p_icol = pattern.icols[i];
3064e1ff37aSBarry Smith     for (j=0; j<p_nnz; j++)  {
307b2f1ab58SBarry Smith       max_row_nnz[i+p_icol[j]]++;
308b2f1ab58SBarry Smith     }
309b2f1ab58SBarry Smith   }
310b2f1ab58SBarry Smith 
311c328eeadSBarry Smith   /* Calculate rows of L */
3124e1ff37aSBarry Smith   for (i = 0; i<nrows; i++)  {
313b2f1ab58SBarry Smith     p_nnz  = pattern.row_nnz[i];
314b2f1ab58SBarry Smith     p_icol = pattern.icols[i];
315b2f1ab58SBarry Smith 
316b2f1ab58SBarry Smith     r_nnz  = retval.row_nnz[i];
317b2f1ab58SBarry Smith     r_icol = retval.icols[i];
318b2f1ab58SBarry Smith     r_val  = retval.values[i];
319b2f1ab58SBarry Smith     A_nnz  = ai[rip[i]+1] - ai[rip[i]];
320b2f1ab58SBarry Smith     A_icol = &aj[ai[rip[i]]];
321b2f1ab58SBarry Smith     A_val  = &aa[ai[rip[i]]];
322b2f1ab58SBarry Smith 
323c328eeadSBarry Smith     /* Calculate  val += A(i,i:n)'; */
3244e1ff37aSBarry Smith     for (j=0; j<A_nnz; j++) {
325b2f1ab58SBarry Smith       k = riip[A_icol[j]];
3262205254eSKarl Rupp       if (k>=i) val[k] = A_val[j];
327b2f1ab58SBarry Smith     }
328b2f1ab58SBarry Smith 
329c328eeadSBarry Smith     /*  Add regularization */
330b2f1ab58SBarry Smith     val[i] += epsdiag;
331b2f1ab58SBarry Smith 
332c328eeadSBarry Smith     /* Calculate lvec   = diag(D(0:i-1)) * L(0:i-1,i);
333c328eeadSBarry Smith         val(i) = A(i,i) - L(0:i-1,i)' * lvec */
3344e1ff37aSBarry Smith     for (j=0; j<r_nnz; j++)  {
335b2f1ab58SBarry Smith       k       = r_icol[j];
336b2f1ab58SBarry Smith       lvec[k] = diag[k] * r_val[j];
337b2f1ab58SBarry Smith       val[i] -= r_val[j] * lvec[k];
338b2f1ab58SBarry Smith     }
339b2f1ab58SBarry Smith 
340c328eeadSBarry Smith     /* Calculate the new diagonal */
341b2f1ab58SBarry Smith     diag[i] = val[i];
3424e1ff37aSBarry Smith     if (PetscRealPart(diag[i])<droptol) {
3439566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL,"Error in spbas_incomplete_cholesky:\n"));
3449566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL,"Negative diagonal in row %" PetscInt_FMT "\n",i+1));
345b2f1ab58SBarry Smith 
346c328eeadSBarry Smith       /* Delete the whole matrix at once. */
3479566063dSJacob Faibussowitsch       PetscCall(spbas_delete(retval));
348cc442ecdSBarry Smith       *success = PETSC_FALSE;
349cc442ecdSBarry Smith       PetscFunctionReturn(0);
350b2f1ab58SBarry Smith     }
351b2f1ab58SBarry Smith 
352c328eeadSBarry Smith     /* If necessary, allocate arrays */
3534e1ff37aSBarry Smith     if (r_nnz==0) {
354cc442ecdSBarry Smith       PetscBool success = spbas_cholesky_row_alloc(retval, i, 1, &n_alloc_used);
35528b400f6SJacob Faibussowitsch       PetscCheck(success,PETSC_COMM_SELF,PETSC_ERR_MEM,"spbas_cholesky_row_alloc() failed");
356b2f1ab58SBarry Smith       r_icol = retval.icols[i];
357b2f1ab58SBarry Smith       r_val  = retval.values[i];
358b2f1ab58SBarry Smith     }
359b2f1ab58SBarry Smith 
360c328eeadSBarry Smith     /* Now, fill in */
361b2f1ab58SBarry Smith     r_icol[r_nnz] = i;
362b2f1ab58SBarry Smith     r_val [r_nnz] = 1.0;
363b2f1ab58SBarry Smith     r_nnz++;
364b2f1ab58SBarry Smith     retval.row_nnz[i]++;
365b2f1ab58SBarry Smith 
366b2f1ab58SBarry Smith     retval.nnz += r_nnz;
367b2f1ab58SBarry Smith 
368c328eeadSBarry Smith     /* Calculate
369c328eeadSBarry Smith         val(i+1:n) = (A(i,i+1:n)- L(0:i-1,i+1:n)' * lvec)/diag(i) */
3704e1ff37aSBarry Smith     for (j=1; j<p_nnz; j++)  {
371b2f1ab58SBarry Smith       k       = i+p_icol[j];
372b2f1ab58SBarry Smith       r1_icol = retval.icols[k];
373b2f1ab58SBarry Smith       r1_val  = retval.values[k];
3744e1ff37aSBarry Smith       for (jL=0; jL<retval.row_nnz[k]; jL++) {
375b2f1ab58SBarry Smith         val[k] -= r1_val[jL] * lvec[r1_icol[jL]];
376b2f1ab58SBarry Smith       }
377b2f1ab58SBarry Smith       val[k] /= diag[i];
378b2f1ab58SBarry Smith 
3794e1ff37aSBarry Smith       if (PetscAbsScalar(val[k]) > droptol || PetscAbsScalar(val[k])< -droptol) {
380c328eeadSBarry Smith         /* If necessary, allocate arrays */
381cc442ecdSBarry Smith         if (!retval.row_nnz[k]) {
382cc442ecdSBarry Smith           PetscBool flag,success = spbas_cholesky_row_alloc(retval, k, max_row_nnz[k], &n_alloc_used);
383cc442ecdSBarry Smith           if (!success) {
3849566063dSJacob Faibussowitsch             PetscCall(spbas_cholesky_garbage_collect(&retval,  i, &n_row_alloc_ok, &n_alloc_used, max_row_nnz));
385cc442ecdSBarry Smith             flag   = spbas_cholesky_row_alloc(retval, k, max_row_nnz[k], &n_alloc_used);
38628b400f6SJacob Faibussowitsch             PetscCheck(flag,PETSC_COMM_SELF,PETSC_ERR_MEM,"Allocation in spbas_cholesky_row_alloc() failed");
387b2f1ab58SBarry Smith             r_icol = retval.icols[i];
388cc442ecdSBarry Smith           }
389b2f1ab58SBarry Smith         }
390b2f1ab58SBarry Smith 
391b2f1ab58SBarry Smith         retval.icols[k][retval.row_nnz[k]]  = i;
392b2f1ab58SBarry Smith         retval.values[k][retval.row_nnz[k]] = val[k];
393b2f1ab58SBarry Smith         retval.row_nnz[k]++;
394b2f1ab58SBarry Smith       }
395b2f1ab58SBarry Smith       val[k] = 0;
396b2f1ab58SBarry Smith     }
397b2f1ab58SBarry Smith 
39871d55d18SBarry Smith     /* Erase the values used in the work arrays */
3992205254eSKarl Rupp     for (j=0; j<r_nnz; j++) lvec[r_icol[j]] = 0;
400b2f1ab58SBarry Smith   }
401b2f1ab58SBarry Smith 
4029566063dSJacob Faibussowitsch   PetscCall(PetscFree(lvec));
4039566063dSJacob Faibussowitsch   PetscCall(PetscFree(val));
404b2f1ab58SBarry Smith 
4059566063dSJacob Faibussowitsch   PetscCall(spbas_cholesky_garbage_collect(&retval, nrows, &n_row_alloc_ok, &n_alloc_used, max_row_nnz));
4069566063dSJacob Faibussowitsch   PetscCall(PetscFree(max_row_nnz));
407b2f1ab58SBarry Smith 
408c328eeadSBarry Smith   /* Place the inverse of the diagonals in the matrix */
4099ccc4a9bSBarry Smith   for (i=0; i<nrows; i++) {
410b2f1ab58SBarry Smith     r_nnz = retval.row_nnz[i];
4112205254eSKarl Rupp 
412b2f1ab58SBarry Smith     retval.values[i][r_nnz-1] = 1.0 / diag[i];
4139ccc4a9bSBarry Smith     for (j=0; j<r_nnz-1; j++) {
414b2f1ab58SBarry Smith       retval.values[i][j] *= -1;
415b2f1ab58SBarry Smith     }
416b2f1ab58SBarry Smith   }
4179566063dSJacob Faibussowitsch   PetscCall(PetscFree(diag));
418b2f1ab58SBarry Smith   *matrix_L = retval;
419cc442ecdSBarry Smith   *success  = PETSC_TRUE;
420b2f1ab58SBarry Smith   PetscFunctionReturn(0);
421b2f1ab58SBarry Smith }
422