xref: /petsc/src/mat/impls/aij/seq/bas/spbas_cholesky.h (revision e32f2f54e699d0aa6e733466c00da7e34666fe5e)
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 */
8b2f1ab58SBarry Smith #undef __FUNCT__
9b2f1ab58SBarry Smith #define __FUNCT__ "spbas_cholesky_row_alloc"
10b2f1ab58SBarry Smith PetscErrorCode spbas_cholesky_row_alloc(spbas_matrix retval, PetscInt k, PetscInt r_nnz,PetscInt * n_alloc_used)
11b2f1ab58SBarry Smith {
12b2f1ab58SBarry Smith    PetscFunctionBegin;
13b2f1ab58SBarry Smith    retval.icols[k]  = &retval.alloc_icol[*n_alloc_used];
14b2f1ab58SBarry Smith    retval.values[k] = &retval.alloc_val[*n_alloc_used];
15b2f1ab58SBarry Smith    *n_alloc_used   += r_nnz;
169ccc4a9bSBarry Smith    if (*n_alloc_used > retval.n_alloc_icol) PetscFunctionReturn(PETSC_ERR_MEM);
17b2f1ab58SBarry Smith    PetscFunctionReturn(0);
18b2f1ab58SBarry Smith }
19b2f1ab58SBarry Smith 
20b2f1ab58SBarry Smith 
21b2f1ab58SBarry Smith /*
22b2f1ab58SBarry Smith   spbas_cholesky_garbage_collect:
23b2f1ab58SBarry Smith      move the rows which have been calculated so far, as well as
24b2f1ab58SBarry Smith      those currently under construction, to the front of the
25b2f1ab58SBarry Smith      array, while putting them in the proper order.
26b2f1ab58SBarry Smith      When it seems necessary, enlarge the current arrays.
27b2f1ab58SBarry Smith 
28b2f1ab58SBarry Smith      NB: re-allocation is being simulated using
29b2f1ab58SBarry Smith          PetscMalloc, memcpy, PetscFree, because
30b2f1ab58SBarry Smith          PetscRealloc does not seem to exist.
31b2f1ab58SBarry Smith 
32b2f1ab58SBarry Smith */
33b2f1ab58SBarry Smith #undef __FUNCT__
34b2f1ab58SBarry Smith #define __FUNCT__ "spbas_cholesky_garbage_collect"
35b2f1ab58SBarry Smith PetscErrorCode spbas_cholesky_garbage_collect(
36c328eeadSBarry Smith       spbas_matrix *result, /* I/O: the Cholesky factor matrix being constructed.
37c328eeadSBarry Smith                                    Only the storage, not the contents of this matrix
38c328eeadSBarry Smith                                    is changed in this function */
39c328eeadSBarry Smith       PetscInt     i_row,             /* I  : Number of rows for which the final contents are
40c328eeadSBarry Smith 				     known */
41c328eeadSBarry Smith       PetscInt     *n_row_alloc_ok,   /* I/O: Number of rows which are already in their final
42c328eeadSBarry Smith 				     places in the arrays: they need not be moved any more */
43c328eeadSBarry Smith       PetscInt     *n_alloc_used,     /* I/O:  */
44c328eeadSBarry Smith       PetscInt     *max_row_nnz       /* I  : Over-estimate of the number of nonzeros needed to
45c328eeadSBarry Smith 				    store each row */
46b2f1ab58SBarry Smith      )
47b2f1ab58SBarry Smith {
48b2f1ab58SBarry Smith 
49b2f1ab58SBarry Smith 
50c328eeadSBarry Smith   /* PSEUDO-CODE:
51c328eeadSBarry Smith     1. Choose the appropriate size for the arrays
52c328eeadSBarry Smith     2. Rescue the arrays which would be lost during garbage collection
53c328eeadSBarry Smith     3. Reallocate and correct administration
54c328eeadSBarry Smith     4. Move all arrays so that they are in proper order */
55b2f1ab58SBarry Smith 
56b2f1ab58SBarry Smith    PetscInt        i,j;
57b2f1ab58SBarry Smith    PetscInt        nrows = result->nrows;
58b2f1ab58SBarry Smith    PetscInt        n_alloc_ok=0;
59b2f1ab58SBarry Smith    PetscInt        n_alloc_ok_max=0;
60b2f1ab58SBarry Smith    PetscErrorCode  ierr;
61b2f1ab58SBarry Smith    PetscInt        need_already= 0;
62b2f1ab58SBarry Smith    PetscInt        n_rows_ahead=0;
63b2f1ab58SBarry Smith    PetscInt        max_need_extra= 0;
64b2f1ab58SBarry Smith    PetscInt        n_alloc_max, n_alloc_est, n_alloc;
65b2f1ab58SBarry Smith    PetscInt        n_alloc_now = result->n_alloc_icol;
66b2f1ab58SBarry Smith    PetscInt        *alloc_icol_old = result->alloc_icol;
67b2f1ab58SBarry Smith    PetscScalar     *alloc_val_old  = result->alloc_val;
68b2f1ab58SBarry Smith    PetscInt        *icol_rescue;
69b2f1ab58SBarry Smith    PetscScalar     *val_rescue;
70b2f1ab58SBarry Smith    PetscInt        n_rescue;
71b2f1ab58SBarry Smith    PetscInt        n_row_rescue;
72b2f1ab58SBarry Smith    PetscInt        i_here, i_last, n_copy;
73d36aa34eSBarry Smith    const PetscReal xtra_perc = 20;
74b2f1ab58SBarry Smith 
75b2f1ab58SBarry Smith    PetscFunctionBegin;
76b2f1ab58SBarry Smith 
77c328eeadSBarry Smith    /*********************************************************
78c328eeadSBarry Smith     1. Choose appropriate array size
79c328eeadSBarry Smith     Count number of rows and memory usage which is already final */
809ccc4a9bSBarry Smith    for (i=0;i<i_row; i++)  {
81b2f1ab58SBarry Smith       n_alloc_ok += result->row_nnz[i];
82b2f1ab58SBarry Smith       n_alloc_ok_max += max_row_nnz[i];
83b2f1ab58SBarry Smith    }
84b2f1ab58SBarry Smith 
85c328eeadSBarry Smith    /* Count currently needed memory usage and future memory requirements
86c328eeadSBarry Smith       (max, predicted)*/
879ccc4a9bSBarry Smith    for (i=i_row; i<nrows; i++) {
889ccc4a9bSBarry Smith      if (!result->row_nnz[i]) {
89b2f1ab58SBarry Smith          max_need_extra  += max_row_nnz[i];
909ccc4a9bSBarry Smith       } else {
91b2f1ab58SBarry Smith          need_already    += max_row_nnz[i];
92b2f1ab58SBarry Smith          n_rows_ahead++;
93b2f1ab58SBarry Smith       }
94b2f1ab58SBarry Smith    }
95b2f1ab58SBarry Smith 
96c328eeadSBarry Smith    /* Make maximal and realistic memory requirement estimates */
97b2f1ab58SBarry Smith    n_alloc_max = n_alloc_ok + need_already + max_need_extra;
989ccc4a9bSBarry Smith    n_alloc_est = n_alloc_ok + need_already + (int) (((PetscReal) max_need_extra) *  ((PetscReal) n_alloc_ok) /((PetscReal) n_alloc_ok_max));
99b2f1ab58SBarry Smith 
100c328eeadSBarry Smith    /* Choose array sizes */
1019ccc4a9bSBarry Smith    if (n_alloc_max == n_alloc_est)  { n_alloc = n_alloc_max; }
1029ccc4a9bSBarry Smith    else if (n_alloc_now >= n_alloc_est)  { n_alloc = n_alloc_now; }
1039ccc4a9bSBarry Smith    else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0))  { n_alloc  = n_alloc_max;  }
1049ccc4a9bSBarry Smith    else { n_alloc = (int) (n_alloc_est * (1+xtra_perc/100.0)); }
105b2f1ab58SBarry Smith 
106c328eeadSBarry Smith    /* If new estimate is less than what we already have,
107c328eeadSBarry Smith       don't reallocate, just garbage-collect */
1089ccc4a9bSBarry Smith    if (n_alloc_max != n_alloc_est && n_alloc < result->n_alloc_icol) {
109b2f1ab58SBarry Smith       n_alloc = result->n_alloc_icol;
110b2f1ab58SBarry Smith    }
111b2f1ab58SBarry Smith 
112c328eeadSBarry Smith    /* Motivate dimension choice */
113c328eeadSBarry Smith    ierr = PetscInfo1(PETSC_NULL,"   Allocating %d nonzeros: ",n_alloc);CHKERRQ(ierr);
1149ccc4a9bSBarry Smith    if (n_alloc_max == n_alloc_est) { ierr = PetscInfo(PETSC_NULL,"this is the correct size\n");CHKERRQ(ierr); }
1159ccc4a9bSBarry Smith    else if (n_alloc_now >= n_alloc_est) { ierr = PetscInfo(PETSC_NULL,"the current size, which seems enough\n");CHKERRQ(ierr); }
1169ccc4a9bSBarry Smith    else if (n_alloc_max < n_alloc_est * (1+xtra_perc/100.0)) { ierr = PetscInfo(PETSC_NULL,"the maximum estimate\n");CHKERRQ(ierr); }
1179ccc4a9bSBarry Smith    else { ierr = PetscInfo1(PETSC_NULL,"%6.2f %% more than the estimate\n",xtra_perc);CHKERRQ(ierr); }
118b2f1ab58SBarry Smith 
119b2f1ab58SBarry Smith 
120c328eeadSBarry Smith    /**********************************************************
121c328eeadSBarry Smith     2. Rescue arrays which would be lost
122c328eeadSBarry Smith     Count how many rows and nonzeros will have to be rescued
123c328eeadSBarry Smith     when moving all arrays in place */
124b2f1ab58SBarry Smith    n_row_rescue = 0; n_rescue = 0;
125b2f1ab58SBarry Smith    if (*n_row_alloc_ok==0) { *n_alloc_used = 0; }
1269ccc4a9bSBarry Smith    else  {
127b2f1ab58SBarry Smith       i = *n_row_alloc_ok - 1;
1289ccc4a9bSBarry Smith       *n_alloc_used = (result->icols[i]-result->alloc_icol) +  result->row_nnz[i];
129b2f1ab58SBarry Smith    }
130b2f1ab58SBarry Smith 
1319ccc4a9bSBarry Smith    for (i=*n_row_alloc_ok; i<nrows; i++) {
132b2f1ab58SBarry Smith       i_here = result->icols[i]-result->alloc_icol;
133b2f1ab58SBarry Smith       i_last = i_here + result->row_nnz[i];
1349ccc4a9bSBarry Smith       if (result->row_nnz[i]>0) {
1359ccc4a9bSBarry Smith          if (*n_alloc_used > i_here || i_last > n_alloc) {
136b2f1ab58SBarry Smith             n_rescue += result->row_nnz[i];
137b2f1ab58SBarry Smith             n_row_rescue++;
138b2f1ab58SBarry Smith          }
139b2f1ab58SBarry Smith 
140b2f1ab58SBarry Smith          if (i<i_row) {  *n_alloc_used += result->row_nnz[i];}
141b2f1ab58SBarry Smith          else         {  *n_alloc_used += max_row_nnz[i];}
142b2f1ab58SBarry Smith       }
143b2f1ab58SBarry Smith    }
144b2f1ab58SBarry Smith 
145c328eeadSBarry Smith    /* Allocate rescue arrays */
1469ccc4a9bSBarry Smith    ierr = PetscMalloc( n_rescue * sizeof(int), &icol_rescue);CHKERRQ(ierr);
1479ccc4a9bSBarry Smith    ierr = PetscMalloc( n_rescue * sizeof(PetscScalar), &val_rescue);CHKERRQ(ierr);
148b2f1ab58SBarry Smith 
149c328eeadSBarry Smith    /* Rescue the arrays which need rescuing */
150b2f1ab58SBarry Smith    n_row_rescue = 0; n_rescue = 0;
151b2f1ab58SBarry Smith    if (*n_row_alloc_ok==0) { *n_alloc_used = 0; }
1529ccc4a9bSBarry Smith    else {
153b2f1ab58SBarry Smith       i = *n_row_alloc_ok - 1;
1549ccc4a9bSBarry Smith       *n_alloc_used = (result->icols[i]-result->alloc_icol) +  result->row_nnz[i];
155b2f1ab58SBarry Smith    }
156b2f1ab58SBarry Smith 
1579ccc4a9bSBarry Smith    for (i=*n_row_alloc_ok; i<nrows; i++) {
158b2f1ab58SBarry Smith       i_here = result->icols[i]-result->alloc_icol;
159b2f1ab58SBarry Smith       i_last = i_here + result->row_nnz[i];
1609ccc4a9bSBarry Smith       if (result->row_nnz[i]>0) {
1614e1ff37aSBarry Smith          if (*n_alloc_used > i_here || i_last > n_alloc) {
1629ccc4a9bSBarry Smith             ierr = PetscMemcpy((void*) &icol_rescue[n_rescue], (void*) result->icols[i], result->row_nnz[i] * sizeof(int));CHKERRQ(ierr);
1639ccc4a9bSBarry Smith             ierr = PetscMemcpy((void*) &val_rescue[n_rescue], (void*) result->values[i], result->row_nnz[i] * sizeof(PetscScalar));CHKERRQ(ierr);
164b2f1ab58SBarry Smith             n_rescue += result->row_nnz[i];
165b2f1ab58SBarry Smith             n_row_rescue++;
166b2f1ab58SBarry Smith          }
167b2f1ab58SBarry Smith 
168b2f1ab58SBarry Smith          if (i<i_row) {  *n_alloc_used += result->row_nnz[i];}
169b2f1ab58SBarry Smith          else         {  *n_alloc_used += max_row_nnz[i];}
170b2f1ab58SBarry Smith       }
171b2f1ab58SBarry Smith    }
172b2f1ab58SBarry Smith 
173c328eeadSBarry Smith    /**********************************************************
174c328eeadSBarry Smith     3. Reallocate and correct administration */
175b2f1ab58SBarry Smith 
1769ccc4a9bSBarry Smith    if (n_alloc != result->n_alloc_icol)  {
177cc6fc633SBarry Smith       n_copy = PetscMin(n_alloc,result->n_alloc_icol);
178b2f1ab58SBarry Smith 
179c328eeadSBarry Smith       /* PETSC knows no REALLOC, so we'll REALLOC ourselves.
180b2f1ab58SBarry Smith 
181c328eeadSBarry Smith 	 Allocate new icol-data, copy old contents */
1829ccc4a9bSBarry Smith       ierr = PetscMalloc(  n_alloc * sizeof(int), &result->alloc_icol);CHKERRQ(ierr);
1839ccc4a9bSBarry Smith       ierr = PetscMemcpy(result->alloc_icol, alloc_icol_old, n_copy*sizeof(int));CHKERRQ(ierr);
184b2f1ab58SBarry Smith 
185c328eeadSBarry Smith       /* Update administration, Reset pointers to new arrays  */
186b2f1ab58SBarry Smith       result->n_alloc_icol = n_alloc;
1879ccc4a9bSBarry Smith       for (i=0; i<nrows; i++) {
1889ccc4a9bSBarry Smith           result->icols[i] =  result->alloc_icol + (result->icols[i]  - alloc_icol_old);
1899ccc4a9bSBarry Smith           result->values[i] =  result->alloc_val + (result->icols[i]  - result->alloc_icol);
190b2f1ab58SBarry Smith       }
191b2f1ab58SBarry Smith 
192c328eeadSBarry Smith       /* Delete old array */
193b2f1ab58SBarry Smith       ierr = PetscFree(alloc_icol_old);CHKERRQ(ierr);
194b2f1ab58SBarry Smith 
195c328eeadSBarry Smith       /* Allocate new value-data, copy old contents */
1969ccc4a9bSBarry Smith       ierr = PetscMalloc(  n_alloc * sizeof(PetscScalar), &result->alloc_val);CHKERRQ(ierr);
1979ccc4a9bSBarry Smith       ierr = PetscMemcpy(result->alloc_val, alloc_val_old, n_copy*sizeof(PetscScalar));CHKERRQ(ierr);
198b2f1ab58SBarry Smith 
199c328eeadSBarry Smith       /* Update administration, Reset pointers to new arrays  */
200b2f1ab58SBarry Smith       result->n_alloc_val  = n_alloc;
2019ccc4a9bSBarry Smith       for (i=0; i<nrows; i++) {
2029ccc4a9bSBarry Smith           result->values[i] =  result->alloc_val + (result->icols[i]  - result->alloc_icol);
203b2f1ab58SBarry Smith       }
204b2f1ab58SBarry Smith 
205c328eeadSBarry Smith       /* Delete old array */
206b2f1ab58SBarry Smith       ierr = PetscFree(alloc_val_old);CHKERRQ(ierr);
207b2f1ab58SBarry Smith    }
208b2f1ab58SBarry Smith 
209b2f1ab58SBarry Smith 
210c328eeadSBarry Smith    /*********************************************************
211c328eeadSBarry Smith     4. Copy all the arrays to their proper places */
212b2f1ab58SBarry Smith    n_row_rescue = 0; n_rescue = 0;
213b2f1ab58SBarry Smith    if (*n_row_alloc_ok==0) { *n_alloc_used = 0; }
2149ccc4a9bSBarry Smith    else {
215b2f1ab58SBarry Smith       i = *n_row_alloc_ok - 1;
2169ccc4a9bSBarry Smith       *n_alloc_used = (result->icols[i]-result->alloc_icol) +  result->row_nnz[i];
217b2f1ab58SBarry Smith    }
218b2f1ab58SBarry Smith 
2199ccc4a9bSBarry Smith    for (i=*n_row_alloc_ok; i<nrows; i++) {
220b2f1ab58SBarry Smith       i_here = result->icols[i]-result->alloc_icol;
221b2f1ab58SBarry Smith       i_last = i_here + result->row_nnz[i];
222b2f1ab58SBarry Smith 
223b2f1ab58SBarry Smith       result->icols[i] = result->alloc_icol + *n_alloc_used;
224b2f1ab58SBarry Smith       result->values[i]= result->alloc_val  + *n_alloc_used;
225b2f1ab58SBarry Smith 
2269ccc4a9bSBarry Smith       if (result->row_nnz[i]>0) {
2279ccc4a9bSBarry Smith          if (*n_alloc_used > i_here || i_last > n_alloc)  {
2289ccc4a9bSBarry Smith             ierr = PetscMemcpy((void*) result->icols[i],  (void*) &icol_rescue[n_rescue], result->row_nnz[i] * sizeof(int));CHKERRQ(ierr);
2299ccc4a9bSBarry Smith             ierr = PetscMemcpy((void*) result->values[i], (void*) &val_rescue[n_rescue],result->row_nnz[i] * sizeof(PetscScalar));CHKERRQ(ierr);
230b2f1ab58SBarry Smith             n_rescue += result->row_nnz[i];
231b2f1ab58SBarry Smith             n_row_rescue++;
2329ccc4a9bSBarry Smith          } else {
2339ccc4a9bSBarry Smith             for (j=0; j<result->row_nnz[i]; j++) {
234b2f1ab58SBarry Smith                 result->icols[i][j]   = result->alloc_icol[i_here+j];
235b2f1ab58SBarry Smith                 result->values[i][j]  = result->alloc_val[ i_here+j];
236b2f1ab58SBarry Smith             }
237b2f1ab58SBarry Smith          }
238b2f1ab58SBarry Smith          if (i<i_row) {  *n_alloc_used += result->row_nnz[i];}
239b2f1ab58SBarry Smith          else         {  *n_alloc_used += max_row_nnz[i];}
240b2f1ab58SBarry Smith       }
241b2f1ab58SBarry Smith    }
242b2f1ab58SBarry Smith 
243c328eeadSBarry Smith    /* Delete the rescue arrays */
244b2f1ab58SBarry Smith    ierr = PetscFree(icol_rescue);CHKERRQ(ierr);
245b2f1ab58SBarry Smith    ierr = PetscFree(val_rescue);CHKERRQ(ierr);
246b2f1ab58SBarry Smith    *n_row_alloc_ok = i_row;
247b2f1ab58SBarry Smith    PetscFunctionReturn(0);
248b2f1ab58SBarry Smith }
249b2f1ab58SBarry Smith 
250b2f1ab58SBarry Smith /*
251b2f1ab58SBarry Smith   spbas_incomplete_cholesky:
252b2f1ab58SBarry Smith      incomplete Cholesky decomposition of a square, symmetric,
253b2f1ab58SBarry Smith      positive definite matrix.
254b2f1ab58SBarry Smith 
255b2f1ab58SBarry Smith      In case negative diagonals are encountered, function returns
256b2f1ab58SBarry Smith      NEGATIVE_DIAGONAL. When this happens, call this function again
257b2f1ab58SBarry Smith      with a larger epsdiag_in, a less sparse pattern, and/or a smaller
258b2f1ab58SBarry Smith      droptol
259b2f1ab58SBarry Smith */
260b2f1ab58SBarry Smith #undef __FUNCT__
261b2f1ab58SBarry Smith #define __FUNCT__ "spbas_incomplete_cholesky"
2629ccc4a9bSBarry 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)
263b2f1ab58SBarry Smith {
264b2f1ab58SBarry Smith    PetscInt         jL;
265b2f1ab58SBarry Smith    Mat_SeqAIJ       *a=(Mat_SeqAIJ*)A->data;
266b2f1ab58SBarry Smith    PetscInt         *ai=a->i,*aj=a->j;
267b2f1ab58SBarry Smith    MatScalar        *aa=a->a;
268b2f1ab58SBarry Smith    PetscInt         nrows, ncols;
269b2f1ab58SBarry Smith    PetscInt         *max_row_nnz;
270b2f1ab58SBarry Smith    PetscErrorCode   ierr;
271b2f1ab58SBarry Smith    spbas_matrix     retval;
272b2f1ab58SBarry Smith    PetscScalar      * diag;
273b2f1ab58SBarry Smith    PetscScalar      * val;
274b2f1ab58SBarry Smith    PetscScalar      * lvec;
275b2f1ab58SBarry Smith    PetscScalar      epsdiag;
276b2f1ab58SBarry Smith    PetscInt         i,j,k;
277b2f1ab58SBarry Smith    const PetscTruth do_values = PETSC_TRUE;
2789ccc4a9bSBarry Smith    PetscInt         * r1_icol;
2799ccc4a9bSBarry Smith    PetscScalar      *r1_val;
2809ccc4a9bSBarry Smith    PetscInt         * r_icol;
2819ccc4a9bSBarry Smith    PetscInt         r_nnz;
2829ccc4a9bSBarry Smith    PetscScalar      *r_val;
2839ccc4a9bSBarry Smith    PetscInt         * A_icol;
2849ccc4a9bSBarry Smith    PetscInt         A_nnz;
2859ccc4a9bSBarry Smith    PetscScalar      *A_val;
2869ccc4a9bSBarry Smith    PetscInt         * p_icol;
2879ccc4a9bSBarry Smith    PetscInt         p_nnz;
2889ccc4a9bSBarry Smith    PetscInt         n_row_alloc_ok = 0;  /* number of rows which have been stored   correctly in the matrix */
2899ccc4a9bSBarry Smith    PetscInt         n_alloc_used = 0;    /* part of result->icols and result->values   which is currently being used */
290b2f1ab58SBarry Smith 
2919ccc4a9bSBarry Smith    PetscFunctionBegin;
2924e1ff37aSBarry Smith    /* Convert the Manteuffel shift from 'fraction of average diagonal' to   dimensioned value */
293b2f1ab58SBarry Smith    ierr = MatGetSize(A, &nrows, &ncols);CHKERRQ(ierr);
294b2f1ab58SBarry Smith    ierr = MatGetTrace(A, &epsdiag);CHKERRQ(ierr);
295b2f1ab58SBarry Smith    epsdiag *= epsdiag_in / nrows;
2969767453cSBarry Smith    ierr = PetscInfo2(PETSC_NULL,"   Dimensioned Manteuffel shift %G Drop tolerance %G\n", PetscRealPart(epsdiag),droptol);CHKERRQ(ierr);
297b2f1ab58SBarry Smith 
298b2f1ab58SBarry Smith    if (droptol<1e-10) {droptol=1e-10;}
299b2f1ab58SBarry Smith 
300*e32f2f54SBarry Smith    if ( (nrows != pattern.nrows) || (ncols != pattern.ncols) || (ncols != nrows) ) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,"Dimension error in spbas_incomplete_cholesky\n");
301b2f1ab58SBarry Smith 
302b2f1ab58SBarry Smith    retval.nrows = nrows;
303b2f1ab58SBarry Smith    retval.ncols = nrows;
304b2f1ab58SBarry Smith    retval.nnz   = pattern.nnz/10;
305b2f1ab58SBarry Smith    retval.col_idx_type = SPBAS_COLUMN_NUMBERS;
306b2f1ab58SBarry Smith    retval.block_data = PETSC_TRUE;
307b2f1ab58SBarry Smith 
308b2f1ab58SBarry Smith    ierr =  spbas_allocate_pattern(&retval, do_values);CHKERRQ(ierr);
3094e1ff37aSBarry Smith    ierr = PetscMemzero((void*) retval.row_nnz, nrows*sizeof(int));CHKERRQ(ierr);
310b2f1ab58SBarry Smith    ierr =  spbas_allocate_data(&retval);CHKERRQ(ierr);
311b2f1ab58SBarry Smith    retval.nnz   = 0;
312b2f1ab58SBarry Smith 
313b2f1ab58SBarry Smith    ierr = PetscMalloc(nrows*sizeof(PetscScalar), &diag);CHKERRQ(ierr);
314b2f1ab58SBarry Smith    ierr = PetscMalloc(nrows*sizeof(PetscScalar), &val);CHKERRQ(ierr);
315b2f1ab58SBarry Smith    ierr = PetscMalloc(nrows*sizeof(PetscScalar), &lvec);CHKERRQ(ierr);
316b2f1ab58SBarry Smith    ierr = PetscMalloc(nrows*sizeof(int), &max_row_nnz);CHKERRQ(ierr);
317*e32f2f54SBarry Smith    if (!(diag && val && lvec && max_row_nnz))  SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MEM, "Allocation error in spbas_incomplete_cholesky\n");
318b2f1ab58SBarry Smith 
3194e1ff37aSBarry Smith    ierr = PetscMemzero((void*) val, nrows*sizeof(PetscScalar));CHKERRQ(ierr);
3204e1ff37aSBarry Smith    ierr = PetscMemzero((void*) lvec, nrows*sizeof(PetscScalar));CHKERRQ(ierr);
3214e1ff37aSBarry Smith    ierr = PetscMemzero((void*) max_row_nnz, nrows*sizeof(int));CHKERRQ(ierr);
322b2f1ab58SBarry Smith 
323c328eeadSBarry Smith    /* Check correct format of sparseness pattern */
324*e32f2f54SBarry Smith    if (pattern.col_idx_type != SPBAS_DIAGONAL_OFFSETS)  SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Error in spbas_incomplete_cholesky: must have diagonal offsets in pattern\n");
325b2f1ab58SBarry Smith 
326c328eeadSBarry Smith    /* Count the nonzeros on transpose of pattern */
3274e1ff37aSBarry Smith    for (i = 0; i<nrows; i++)  {
328b2f1ab58SBarry Smith       p_nnz  = pattern.row_nnz[i];
329b2f1ab58SBarry Smith       p_icol = pattern.icols[i];
3304e1ff37aSBarry Smith       for (j=0; j<p_nnz; j++)  {
331b2f1ab58SBarry Smith          max_row_nnz[i+p_icol[j]]++;
332b2f1ab58SBarry Smith       }
333b2f1ab58SBarry Smith    }
334b2f1ab58SBarry Smith 
335c328eeadSBarry Smith    /* Calculate rows of L */
3364e1ff37aSBarry Smith    for (i = 0; i<nrows; i++)  {
337b2f1ab58SBarry Smith       p_nnz  = pattern.row_nnz[i];
338b2f1ab58SBarry Smith       p_icol = pattern.icols[i];
339b2f1ab58SBarry Smith 
340b2f1ab58SBarry Smith       r_nnz  = retval.row_nnz[i];
341b2f1ab58SBarry Smith       r_icol = retval.icols[i];
342b2f1ab58SBarry Smith       r_val  = retval.values[i];
343b2f1ab58SBarry Smith       A_nnz  = ai[rip[i]+1] - ai[rip[i]];
344b2f1ab58SBarry Smith       A_icol = &aj[ai[rip[i]]];
345b2f1ab58SBarry Smith       A_val  = &aa[ai[rip[i]]];
346b2f1ab58SBarry Smith 
347c328eeadSBarry Smith       /* Calculate  val += A(i,i:n)'; */
3484e1ff37aSBarry Smith       for (j=0; j<A_nnz; j++) {
349b2f1ab58SBarry Smith          k = riip[A_icol[j]];
350b2f1ab58SBarry Smith          if (k>=i) { val[k] = A_val[j]; }
351b2f1ab58SBarry Smith       }
352b2f1ab58SBarry Smith 
353c328eeadSBarry Smith       /*  Add regularization */
354b2f1ab58SBarry Smith       val[i] += epsdiag;
355b2f1ab58SBarry Smith 
356c328eeadSBarry Smith       /* Calculate lvec   = diag(D(0:i-1)) * L(0:i-1,i);
357c328eeadSBarry Smith 	 val(i) = A(i,i) - L(0:i-1,i)' * lvec */
3584e1ff37aSBarry Smith       for ( j=0; j<r_nnz; j++)  {
359b2f1ab58SBarry Smith          k           = r_icol[j];
360b2f1ab58SBarry Smith          lvec[k]     = diag[k] * r_val[j];
361b2f1ab58SBarry Smith          val[i]     -= r_val[j] * lvec[k];
362b2f1ab58SBarry Smith       }
363b2f1ab58SBarry Smith 
364c328eeadSBarry Smith       /* Calculate the new diagonal */
365b2f1ab58SBarry Smith       diag[i] = val[i];
3664e1ff37aSBarry Smith       if (PetscRealPart(diag[i])<droptol)  {
367c328eeadSBarry Smith          ierr = PetscInfo(PETSC_NULL,"Error in spbas_incomplete_cholesky:\n");
368c328eeadSBarry Smith          ierr = PetscInfo1(PETSC_NULL,"Negative diagonal in row %d\n",i+1);
369b2f1ab58SBarry Smith 
370c328eeadSBarry Smith          /* Delete the whole matrix at once. */
371b2f1ab58SBarry Smith          ierr = spbas_delete(retval);
372b2f1ab58SBarry Smith          return NEGATIVE_DIAGONAL;
373b2f1ab58SBarry Smith       }
374b2f1ab58SBarry Smith 
375c328eeadSBarry Smith       /* If necessary, allocate arrays */
3764e1ff37aSBarry Smith       if (r_nnz==0) {
3779ccc4a9bSBarry Smith          ierr = spbas_cholesky_row_alloc( retval, i, 1, &n_alloc_used);CHKERRQ(ierr);
3784e1ff37aSBarry Smith          if (ierr == PETSC_ERR_MEM) {
3799ccc4a9bSBarry Smith             ierr = spbas_cholesky_garbage_collect( &retval,  i, &n_row_alloc_ok, &n_alloc_used, max_row_nnz);CHKERRQ(ierr);
380b2f1ab58SBarry Smith             ierr = spbas_cholesky_row_alloc( retval, i, 1, &n_alloc_used);
381b2f1ab58SBarry Smith          }
382b2f1ab58SBarry Smith          r_icol = retval.icols[i];
383b2f1ab58SBarry Smith          r_val  = retval.values[i];
384b2f1ab58SBarry Smith       }
385b2f1ab58SBarry Smith 
386c328eeadSBarry Smith       /* Now, fill in */
387b2f1ab58SBarry Smith       r_icol[r_nnz] = i;
388b2f1ab58SBarry Smith       r_val [r_nnz] = 1.0;
389b2f1ab58SBarry Smith       r_nnz++;
390b2f1ab58SBarry Smith       retval.row_nnz[i]++;
391b2f1ab58SBarry Smith 
392b2f1ab58SBarry Smith       retval.nnz += r_nnz;
393b2f1ab58SBarry Smith 
394c328eeadSBarry Smith       /* Calculate
395c328eeadSBarry Smith          val(i+1:n) = (A(i,i+1:n)- L(0:i-1,i+1:n)' * lvec)/diag(i) */
3964e1ff37aSBarry Smith       for (j=1; j<p_nnz; j++)  {
397b2f1ab58SBarry Smith          k = i+p_icol[j];
398b2f1ab58SBarry Smith          r1_icol = retval.icols[k];
399b2f1ab58SBarry Smith          r1_val  = retval.values[k];
4004e1ff37aSBarry Smith          for (jL=0; jL<retval.row_nnz[k]; jL++)   {
401b2f1ab58SBarry Smith             val[k] -= r1_val[jL] * lvec[r1_icol[jL]];
402b2f1ab58SBarry Smith          }
403b2f1ab58SBarry Smith          val[k] /= diag[i];
404b2f1ab58SBarry Smith 
4054e1ff37aSBarry Smith          if (PetscAbsScalar(val[k]) > droptol || PetscAbsScalar(val[k])< -droptol) {
406c328eeadSBarry Smith 	   /* If necessary, allocate arrays */
4074e1ff37aSBarry Smith             if (retval.row_nnz[k]==0) {
4089ccc4a9bSBarry Smith                ierr = spbas_cholesky_row_alloc( retval, k, max_row_nnz[k], &n_alloc_used);
4099ccc4a9bSBarry Smith                if (ierr == PETSC_ERR_MEM) {
4109ccc4a9bSBarry Smith                   ierr = spbas_cholesky_garbage_collect( &retval,  i, &n_row_alloc_ok, &n_alloc_used, max_row_nnz);CHKERRQ(ierr);
4119ccc4a9bSBarry Smith                   ierr = spbas_cholesky_row_alloc( retval, k, max_row_nnz[k], &n_alloc_used);CHKERRQ(ierr);
412b2f1ab58SBarry Smith                   r_icol = retval.icols[i];
413b2f1ab58SBarry Smith                   r_val  = retval.values[i];
414b2f1ab58SBarry Smith                }
415b2f1ab58SBarry Smith             }
416b2f1ab58SBarry Smith 
417b2f1ab58SBarry Smith             retval.icols[k][retval.row_nnz[k]] = i;
418b2f1ab58SBarry Smith             retval.values[k][retval.row_nnz[k]] = val[k];
419b2f1ab58SBarry Smith             retval.row_nnz[k]++;
420b2f1ab58SBarry Smith          }
421b2f1ab58SBarry Smith          val[k] = 0;
422b2f1ab58SBarry Smith       }
423b2f1ab58SBarry Smith 
42471d55d18SBarry Smith       /* Erase the values used in the work arrays */
425b2f1ab58SBarry Smith       for (j=0; j<r_nnz; j++) { lvec[r_icol[j]] = 0; }
426b2f1ab58SBarry Smith    }
427b2f1ab58SBarry Smith 
4289ccc4a9bSBarry Smith    ierr=PetscFree(lvec);CHKERRQ(ierr);
4299ccc4a9bSBarry Smith    ierr=PetscFree(val);CHKERRQ(ierr);
430b2f1ab58SBarry Smith 
4319ccc4a9bSBarry Smith    ierr = spbas_cholesky_garbage_collect( &retval, nrows, &n_row_alloc_ok, &n_alloc_used, max_row_nnz);CHKERRQ(ierr);
432b2f1ab58SBarry Smith    ierr=PetscFree(max_row_nnz);CHKERRQ(ierr);
433b2f1ab58SBarry Smith 
434c328eeadSBarry Smith    /* Place the inverse of the diagonals in the matrix */
4359ccc4a9bSBarry Smith    for (i=0; i<nrows; i++) {
436b2f1ab58SBarry Smith       r_nnz = retval.row_nnz[i];
437b2f1ab58SBarry Smith       retval.values[i][r_nnz-1] = 1.0 / diag[i];
4389ccc4a9bSBarry Smith       for (j=0; j<r_nnz-1; j++) {
439b2f1ab58SBarry Smith          retval.values[i][j] *= -1;
440b2f1ab58SBarry Smith       }
441b2f1ab58SBarry Smith    }
442b2f1ab58SBarry Smith    ierr=PetscFree(diag);CHKERRQ(ierr);
443b2f1ab58SBarry Smith    *matrix_L = retval;
444b2f1ab58SBarry Smith    PetscFunctionReturn(0);
445b2f1ab58SBarry Smith }
446b2f1ab58SBarry Smith 
447