xref: /petsc/src/mat/impls/aij/seq/bas/spbas.c (revision 0ad7597df3318518cf08ebbaa20b6a237bca446d)
1 #include <../src/mat/impls/aij/seq/aij.h>
2 #include <../src/mat/impls/aij/seq/bas/spbas.h>
3 
4 /*MC
5   MATSOLVERBAS -  Provides ICC(k) with drop tolerance
6 
7   Works with MATAIJ  matrices
8 
9   Options Database Keys:
10 + -pc_factor_levels <l>
11 - -pc_factor_drop_tolerance
12 
13   Level: beginner
14 
15    Contributed by: Bas van 't Hof
16 
17 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, PCFactorSetLevels(), PCFactorSetDropTolerance()
18 
19 M*/
20 
21 /*
22   spbas_memory_requirement:
23     Calculate the number of bytes needed to store tha matrix
24 */
25 #undef __FUNCT__
26 #define __FUNCT__ "spbas_memory_requirement"
27 long int spbas_memory_requirement(spbas_matrix matrix)
28 {
29    long int memreq = 6 * sizeof(PetscInt)         + /* nrows, ncols, nnz, n_alloc_icol,
30                                                        n_alloc_val, col_idx_type */
31      sizeof(PetscBool)                 + /* block_data */
32      sizeof(PetscScalar**)        + /* values */
33      sizeof(PetscScalar*)         + /* alloc_val */
34      2 * sizeof(int**)            + /* icols, icols0 */
35      2 * sizeof(PetscInt*)             + /* row_nnz, alloc_icol */
36      matrix.nrows * sizeof(PetscInt)   + /* row_nnz[*] */
37      matrix.nrows * sizeof(PetscInt*);   /* icols[*] */
38 
39    /* icol0[*] */
40    if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) { memreq += matrix.nrows * sizeof(PetscInt); }
41 
42    /* icols[*][*] */
43    if (matrix.block_data) { memreq += matrix.n_alloc_icol * sizeof(PetscInt); }
44    else { memreq += matrix.nnz * sizeof(PetscInt); }
45 
46    if (matrix.values) {
47      memreq += matrix.nrows * sizeof(PetscScalar*); /* values[*] */
48      /* values[*][*] */
49        if (matrix.block_data) { memreq += matrix.n_alloc_val  * sizeof(PetscScalar); }
50        else { memreq += matrix.nnz * sizeof(PetscScalar); }
51    }
52    return memreq;
53 }
54 
55 /*
56   spbas_allocate_pattern:
57     allocate the pattern arrays row_nnz, icols and optionally values
58 */
59 #undef __FUNCT__
60 #define __FUNCT__ "spbas_allocate_pattern"
61 PetscErrorCode spbas_allocate_pattern(spbas_matrix * result, PetscBool  do_values)
62 {
63    PetscErrorCode ierr;
64    PetscInt       nrows = result->nrows;
65    PetscInt       col_idx_type = result->col_idx_type;
66 
67    PetscFunctionBegin;
68    /* Allocate sparseness pattern */
69    ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->row_nnz);CHKERRQ(ierr);
70    ierr = PetscMalloc(nrows*sizeof(PetscInt*),&result->icols);CHKERRQ(ierr);
71 
72    /* If offsets are given wrt an array, create array */
73    if (col_idx_type == SPBAS_OFFSET_ARRAY) {
74       ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->icol0);CHKERRQ(ierr);
75    } else  {
76       result->icol0 = PETSC_NULL;
77    }
78 
79    /* If values are given, allocate values array */
80    if (do_values)  {
81       ierr = PetscMalloc(nrows*sizeof(PetscScalar*),&result->values);CHKERRQ(ierr);
82    } else {
83       result->values = PETSC_NULL;
84    }
85    PetscFunctionReturn(0);
86 }
87 
88 /*
89 spbas_allocate_data:
90    in case of block_data:
91        Allocate the data arrays alloc_icol and optionally alloc_val,
92        set appropriate pointers from icols and values;
93    in case of !block_data:
94        Allocate the arrays icols[i] and optionally values[i]
95 */
96 #undef __FUNCT__
97 #define __FUNCT__ "spbas_allocate_data"
98 PetscErrorCode spbas_allocate_data(spbas_matrix * result)
99 {
100    PetscInt       i;
101    PetscInt       nnz   = result->nnz;
102    PetscInt       nrows = result->nrows;
103    PetscInt       r_nnz;
104    PetscErrorCode ierr;
105    PetscBool      do_values = (result->values != PETSC_NULL) ? PETSC_TRUE : PETSC_FALSE;
106    PetscBool      block_data = result->block_data;
107 
108    PetscFunctionBegin;
109    if (block_data) {
110      /* Allocate the column number array and point to it */
111       result->n_alloc_icol = nnz;
112       ierr = PetscMalloc(nnz*sizeof(PetscInt), &result->alloc_icol);CHKERRQ(ierr);
113       result->icols[0]= result->alloc_icol;
114       for (i=1; i<nrows; i++)  {
115           result->icols[i] = result->icols[i-1] + result->row_nnz[i-1];
116       }
117 
118       /* Allocate the value array and point to it */
119       if (do_values)  {
120          result->n_alloc_val= nnz;
121          ierr = PetscMalloc(nnz*sizeof(PetscScalar), &result->alloc_val);CHKERRQ(ierr);
122          result->values[0]= result->alloc_val;
123          for (i=1; i<nrows; i++) {
124              result->values[i] = result->values[i-1] + result->row_nnz[i-1];
125          }
126       }
127    } else {
128       for (i=0; i<nrows; i++)   {
129          r_nnz = result->row_nnz[i];
130          ierr = PetscMalloc(r_nnz*sizeof(PetscInt), &result->icols[i]);CHKERRQ(ierr);
131       }
132       if (do_values) {
133          for (i=0; i<nrows; i++)  {
134             r_nnz = result->row_nnz[i];
135             ierr = PetscMalloc(r_nnz*sizeof(PetscScalar), &result->values[i]);CHKERRQ(ierr);
136          }
137       }
138    }
139    PetscFunctionReturn(0);
140 }
141 
142 /*
143   spbas_row_order_icol
144      determine if row i1 should come
145        + before row i2 in the sorted rows (return -1),
146        + after (return 1)
147        + is identical (return 0).
148 */
149 #undef __FUNCT__
150 #define __FUNCT__ "spbas_row_order_icol"
151 int spbas_row_order_icol(PetscInt i1, PetscInt i2, PetscInt *irow_in, PetscInt *icol_in,PetscInt col_idx_type)
152 {
153    PetscInt j;
154    PetscInt nnz1 = irow_in[i1+1] - irow_in[i1];
155    PetscInt nnz2 = irow_in[i2+1] - irow_in[i2];
156    PetscInt * icol1 = &icol_in[irow_in[i1]];
157    PetscInt * icol2 = &icol_in[irow_in[i2]];
158 
159    if (nnz1<nnz2) {return -1;}
160    if (nnz1>nnz2) {return 1;}
161 
162    if (col_idx_type == SPBAS_COLUMN_NUMBERS) {
163       for (j=0; j<nnz1; j++) {
164          if (icol1[j]< icol2[j]) {return -1;}
165          if (icol1[j]> icol2[j]) {return 1;}
166       }
167    } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {
168       for (j=0; j<nnz1; j++) {
169          if (icol1[j]-i1< icol2[j]-i2) {return -1;}
170          if (icol1[j]-i1> icol2[j]-i2) {return 1;}
171       }
172    } else if (col_idx_type == SPBAS_OFFSET_ARRAY) {
173       for (j=1; j<nnz1; j++) {
174          if (icol1[j]-icol1[0] < icol2[j]-icol2[0]) {return -1;}
175          if (icol1[j]-icol1[0] > icol2[j]-icol2[0]) {return 1;}
176       }
177    }
178    return 0;
179 }
180 
181 /*
182   spbas_mergesort_icols:
183     return a sorting of the rows in which identical sparseness patterns are
184     next to each other
185 */
186 #undef __FUNCT__
187 #define __FUNCT__ "spbas_mergesort_icols"
188 PetscErrorCode spbas_mergesort_icols(PetscInt nrows, PetscInt * irow_in, PetscInt * icol_in,PetscInt col_idx_type, PetscInt *isort)
189 {
190    PetscErrorCode ierr;
191    PetscInt       istep;       /* Chunk-sizes of already sorted parts of arrays */
192    PetscInt       i, i1, i2;   /* Loop counters for (partly) sorted arrays */
193    PetscInt       istart, i1end, i2end; /* start of newly sorted array part, end of both  parts */
194    PetscInt       *ialloc;     /* Allocated arrays */
195    PetscInt       *iswap;      /* auxiliary pointers for swapping */
196    PetscInt       *ihlp1;      /* Pointers to new version of arrays, */
197    PetscInt       *ihlp2;      /* Pointers to previous version of arrays, */
198 
199    PetscFunctionBegin;
200    ierr = PetscMalloc(nrows*sizeof(PetscInt),&ialloc);CHKERRQ(ierr);
201 
202    ihlp1 = ialloc;
203    ihlp2 = isort;
204 
205    /* Sorted array chunks are first 1 long, and increase until they are the complete array */
206    for (istep=1; istep<nrows; istep*=2) {
207      /*
208        Combine sorted parts
209             istart:istart+istep-1 and istart+istep-1:istart+2*istep-1
210        of ihlp2 and vhlp2
211 
212        into one sorted part
213             istart:istart+2*istep-1
214        of ihlp1 and vhlp1
215      */
216       for (istart=0; istart<nrows; istart+=2*istep) {
217 
218          /* Set counters and bound array part endings */
219          i1=istart;        i1end = i1+istep;  if (i1end>nrows) {i1end=nrows;}
220          i2=istart+istep;  i2end = i2+istep;  if (i2end>nrows) {i2end=nrows;}
221 
222          /* Merge the two array parts */
223          for (i=istart; i<i2end; i++)  {
224             if (i1<i1end && i2<i2end && spbas_row_order_icol(ihlp2[i1], ihlp2[i2], irow_in, icol_in, col_idx_type) < 0) {
225                 ihlp1[i] = ihlp2[i1];
226                 i1++;
227             } else if (i2<i2end) {
228                 ihlp1[i] = ihlp2[i2];
229                 i2++;
230             } else {
231                 ihlp1[i] = ihlp2[i1];
232                 i1++;
233             }
234          }
235       }
236 
237       /* Swap the two array sets */
238       iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap;
239    }
240 
241    /* Copy one more time in case the sorted arrays are the temporary ones */
242    if (ihlp2 != isort) {
243       for (i=0; i<nrows; i++) { isort[i] = ihlp2[i]; }
244    }
245    ierr = PetscFree(ialloc);CHKERRQ(ierr);
246    PetscFunctionReturn(0);
247 }
248 
249 
250 
251 /*
252   spbas_compress_pattern:
253      calculate a compressed sparseness pattern for a sparseness pattern
254      given in compressed row storage. The compressed sparseness pattern may
255      require (much) less memory.
256 */
257 #undef __FUNCT__
258 #define __FUNCT__ "spbas_compress_pattern"
259 PetscErrorCode spbas_compress_pattern(PetscInt *irow_in, PetscInt *icol_in, PetscInt nrows, PetscInt ncols, PetscInt col_idx_type, spbas_matrix *B,PetscReal *mem_reduction)
260 {
261    PetscInt         nnz = irow_in[nrows];
262    long int         mem_orig = (nrows + nnz) * sizeof(PetscInt);
263    long int         mem_compressed;
264    PetscErrorCode   ierr;
265    PetscInt         *isort;
266    PetscInt         *icols;
267    PetscInt         row_nnz;
268    PetscInt         *ipoint;
269    PetscBool        *used;
270    PetscInt         ptr;
271    PetscInt         i,j;
272    const PetscBool  no_values = PETSC_FALSE;
273 
274    PetscFunctionBegin;
275    /* Allocate the structure of the new matrix */
276    B->nrows = nrows;
277    B->ncols = ncols;
278    B->nnz   = nnz;
279    B->col_idx_type= col_idx_type;
280    B->block_data = PETSC_TRUE;
281    ierr = spbas_allocate_pattern(B, no_values);CHKERRQ(ierr);
282 
283    /* When using an offset array, set it */
284    if (col_idx_type==SPBAS_OFFSET_ARRAY)  {
285       for (i=0; i<nrows; i++) {B->icol0[i] = icol_in[irow_in[i]];}
286    }
287 
288    /* Allocate the ordering for the rows */
289    ierr = PetscMalloc(nrows*sizeof(PetscInt),&isort);CHKERRQ(ierr);
290    ierr = PetscMalloc(nrows*sizeof(PetscInt),&ipoint);CHKERRQ(ierr);
291    ierr = PetscMalloc(nrows*sizeof(PetscBool),&used);CHKERRQ(ierr);
292 
293    /*  Initialize the sorting */
294    ierr = PetscMemzero((void*) used, nrows*sizeof(PetscBool));CHKERRQ(ierr);
295    for (i = 0; i<nrows; i++)  {
296       B->row_nnz[i] = irow_in[i+1]-irow_in[i];
297       isort[i] = i;
298       ipoint[i]= i;
299    }
300 
301    /* Sort the rows so that identical columns will be next to each other */
302    ierr = spbas_mergesort_icols(nrows, irow_in, icol_in, col_idx_type, isort);CHKERRQ(ierr);
303    ierr = PetscInfo(PETSC_NULL,"Rows have been sorted for patterns\n");CHKERRQ(ierr);
304 
305    /* Replace identical rows with the first one in the list */
306    for (i=1; i<nrows; i++) {
307       if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) {
308          ipoint[isort[i]] = ipoint[isort[i-1]];
309       }
310    }
311 
312    /* Collect the rows which are used*/
313    for (i=0; i<nrows; i++) {used[ipoint[i]] = PETSC_TRUE;}
314 
315    /* Calculate needed memory */
316    B->n_alloc_icol = 0;
317    for (i=0; i<nrows; i++)  {
318       if (used[i]) {B->n_alloc_icol += B->row_nnz[i];}
319    }
320    ierr = PetscMalloc(B->n_alloc_icol*sizeof(PetscInt),&B->alloc_icol);CHKERRQ(ierr);
321 
322    /* Fill in the diagonal offsets for the rows which store their own data */
323    ptr = 0;
324    for (i=0; i<B->nrows; i++) {
325       if (used[i])  {
326          B->icols[i] = &B->alloc_icol[ptr];
327          icols = &icol_in[irow_in[i]];
328          row_nnz = B->row_nnz[i];
329          if (col_idx_type == SPBAS_COLUMN_NUMBERS) {
330             for (j=0; j<row_nnz; j++) {
331                B->icols[i][j] = icols[j];
332             }
333          } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {
334             for (j=0; j<row_nnz; j++) {
335                B->icols[i][j] = icols[j]-i;
336             }
337          } else if (col_idx_type == SPBAS_OFFSET_ARRAY) {
338             for (j=0; j<row_nnz; j++) {
339                B->icols[i][j] = icols[j]-icols[0];
340             }
341          }
342          ptr += B->row_nnz[i];
343       }
344    }
345 
346    /* Point to the right places for all data */
347    for (i=0; i<nrows; i++) {
348       B->icols[i] = B->icols[ipoint[i]];
349    }
350    ierr = PetscInfo(PETSC_NULL,"Row patterns have been compressed\n");CHKERRQ(ierr);
351    ierr = PetscInfo1(PETSC_NULL,"         (%G nonzeros per row)\n",  (PetscReal) nnz / (PetscReal) nrows);CHKERRQ(ierr);
352 
353    ierr=PetscFree(isort);CHKERRQ(ierr);
354    ierr=PetscFree(used);CHKERRQ(ierr);
355    ierr=PetscFree(ipoint);CHKERRQ(ierr);
356 
357    mem_compressed = spbas_memory_requirement(*B);
358    *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig;
359    PetscFunctionReturn(0);
360 }
361 
362 /*
363    spbas_incomplete_cholesky
364        Incomplete Cholesky decomposition
365 */
366 #include <../src/mat/impls/aij/seq/bas/spbas_cholesky.h>
367 
368 /*
369   spbas_delete : de-allocate the arrays owned by this matrix
370 */
371 #undef __FUNCT__
372 #define __FUNCT__ "spbas_delete"
373 PetscErrorCode spbas_delete(spbas_matrix matrix)
374 {
375    PetscInt       i;
376    PetscErrorCode ierr;
377 
378    PetscFunctionBegin;
379    if (matrix.block_data) {
380       ierr=PetscFree(matrix.alloc_icol);CHKERRQ(ierr);
381       if (matrix.values) {ierr=PetscFree(matrix.alloc_val);CHKERRQ(ierr);}
382    } else  {
383      for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.icols[i]);CHKERRQ(ierr);}
384      ierr = PetscFree(matrix.icols);CHKERRQ(ierr);
385      if (matrix.values) {
386          for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.values[i]);CHKERRQ(ierr);}
387       }
388    }
389 
390    ierr=PetscFree(matrix.row_nnz);CHKERRQ(ierr);
391    ierr=PetscFree(matrix.icols);CHKERRQ(ierr);
392    if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) {ierr=PetscFree(matrix.icol0);CHKERRQ(ierr);}
393    ierr=PetscFree(matrix.values);CHKERRQ(ierr);
394    PetscFunctionReturn(0);
395 }
396 
397 /*
398 spbas_matrix_to_crs:
399    Convert an spbas_matrix to compessed row storage
400 */
401 #undef __FUNCT__
402 #define __FUNCT__ "spbas_matrix_to_crs"
403 PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out)
404 {
405    PetscInt       nrows = matrix_A.nrows;
406    PetscInt       nnz   = matrix_A.nnz;
407    PetscInt       i,j,r_nnz,i0;
408    PetscInt       *irow;
409    PetscInt       *icol;
410    PetscInt       *icol_A;
411    MatScalar      *val;
412    PetscScalar    *val_A;
413    PetscInt       col_idx_type = matrix_A.col_idx_type;
414    PetscBool      do_values = matrix_A.values ? PETSC_TRUE : PETSC_FALSE;
415    PetscErrorCode ierr;
416 
417    PetscFunctionBegin;
418    ierr = PetscMalloc(sizeof(PetscInt) * (nrows+1), &irow);CHKERRQ(ierr);
419    ierr = PetscMalloc(sizeof(PetscInt) * nnz, &icol);CHKERRQ(ierr);
420    *icol_out = icol;
421    *irow_out=irow;
422    if (do_values)  {
423       ierr = PetscMalloc(sizeof(MatScalar) * nnz, &val);CHKERRQ(ierr);
424       *val_out = val; *icol_out = icol; *irow_out=irow;
425    }
426 
427    irow[0]=0;
428    for (i=0; i<nrows; i++) {
429       r_nnz = matrix_A.row_nnz[i];
430       i0 = irow[i];
431       irow[i+1] = i0 + r_nnz;
432       icol_A = matrix_A.icols[i];
433 
434       if (do_values) {
435           val_A = matrix_A.values[i];
436           for (j=0; j<r_nnz; j++) {
437              icol[i0+j] = icol_A[j];
438              val[i0+j]  = val_A[j];
439           }
440       } else {
441           for (j=0; j<r_nnz; j++) { icol[i0+j] = icol_A[j]; }
442       }
443 
444       if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {
445          for (j=0; j<r_nnz; j++) { icol[i0+j] += i; }
446       }
447       else if (col_idx_type == SPBAS_OFFSET_ARRAY) {
448          i0 = matrix_A.icol0[i];
449          for (j=0; j<r_nnz; j++) { icol[i0+j] += i0; }
450       }
451    }
452    PetscFunctionReturn(0);
453 }
454 
455 
456 /*
457     spbas_transpose
458        return the transpose of a matrix
459 */
460 #undef __FUNCT__
461 #define __FUNCT__ "spbas_transpose"
462 PetscErrorCode spbas_transpose(spbas_matrix in_matrix, spbas_matrix * result)
463 {
464    PetscInt       col_idx_type = in_matrix.col_idx_type;
465    PetscInt       nnz   = in_matrix.nnz;
466    PetscInt       ncols = in_matrix.nrows;
467    PetscInt       nrows = in_matrix.ncols;
468    PetscInt       i,j,k;
469    PetscInt       r_nnz;
470    PetscInt       *irow;
471    PetscInt       icol0 = 0;
472    PetscScalar    * val;
473    PetscErrorCode ierr;
474 
475    PetscFunctionBegin;
476    /* Copy input values */
477    result->nrows        = nrows;
478    result->ncols        = ncols;
479    result->nnz          = nnz;
480    result->col_idx_type = SPBAS_COLUMN_NUMBERS;
481    result->block_data   = PETSC_TRUE;
482 
483    /* Allocate sparseness pattern */
484    ierr =  spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
485 
486    /*  Count the number of nonzeros in each row */
487    for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; }
488 
489    for (i=0; i<ncols; i++) {
490       r_nnz = in_matrix.row_nnz[i];
491       irow  = in_matrix.icols[i];
492       if (col_idx_type == SPBAS_COLUMN_NUMBERS)  {
493          for (j=0; j<r_nnz; j++) { result->row_nnz[irow[j]]++; }
494       } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS)  {
495          for (j=0; j<r_nnz; j++) { result->row_nnz[i+irow[j]]++; }
496       } else if (col_idx_type == SPBAS_OFFSET_ARRAY) {
497          icol0=in_matrix.icol0[i];
498          for (j=0; j<r_nnz; j++) { result->row_nnz[icol0+irow[j]]++; }
499       }
500    }
501 
502    /* Set the pointers to the data */
503    ierr = spbas_allocate_data(result);CHKERRQ(ierr);
504 
505    /* Reset the number of nonzeros in each row */
506    for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; }
507 
508    /* Fill the data arrays */
509    if (in_matrix.values) {
510       for (i=0; i<ncols; i++) {
511          r_nnz = in_matrix.row_nnz[i];
512          irow  = in_matrix.icols[i];
513          val   = in_matrix.values[i];
514 
515          if      (col_idx_type == SPBAS_COLUMN_NUMBERS)   {icol0=0;}
516          else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;}
517          else if (col_idx_type == SPBAS_OFFSET_ARRAY)     {icol0=in_matrix.icol0[i];}
518          for (j=0; j<r_nnz; j++)  {
519             k = icol0 + irow[j];
520             result->icols[k][result->row_nnz[k]]  = i;
521             result->values[k][result->row_nnz[k]] = val[j];
522             result->row_nnz[k]++;
523          }
524       }
525    }  else {
526       for (i=0; i<ncols; i++) {
527          r_nnz = in_matrix.row_nnz[i];
528          irow  = in_matrix.icols[i];
529 
530          if      (col_idx_type == SPBAS_COLUMN_NUMBERS)   {icol0=0;}
531          else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;}
532          else if (col_idx_type == SPBAS_OFFSET_ARRAY)     {icol0=in_matrix.icol0[i];}
533 
534          for (j=0; j<r_nnz; j++) {
535             k = icol0 + irow[j];
536             result->icols[k][result->row_nnz[k]]  = i;
537             result->row_nnz[k]++;
538          }
539       }
540    }
541    PetscFunctionReturn(0);
542 }
543 
544 /*
545    spbas_mergesort
546 
547       mergesort for an array of intergers and an array of associated
548       reals
549 
550       on output, icol[0..nnz-1] is increasing;
551                   val[0..nnz-1] has undergone the same permutation as icol
552 
553       NB: val may be PETSC_NULL: in that case, only the integers are sorted
554 
555 */
556 #undef __FUNCT__
557 #define __FUNCT__ "spbas_mergesort"
558 PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val)
559 {
560   PetscInt       istep;       /* Chunk-sizes of already sorted parts of arrays */
561   PetscInt       i, i1, i2;   /* Loop counters for (partly) sorted arrays */
562   PetscInt       istart, i1end, i2end; /* start of newly sorted array part, end of both parts */
563   PetscInt       *ialloc;     /* Allocated arrays */
564   PetscScalar    *valloc=PETSC_NULL;
565   PetscInt       *iswap;      /* auxiliary pointers for swapping */
566   PetscScalar    *vswap;
567   PetscInt       *ihlp1;      /* Pointers to new version of arrays, */
568   PetscScalar    *vhlp1=PETSC_NULL;  /* (arrays under construction) */
569   PetscInt       *ihlp2;      /* Pointers to previous version of arrays, */
570   PetscScalar    *vhlp2=PETSC_NULL;
571   PetscErrorCode ierr;
572 
573    ierr = PetscMalloc(nnz*sizeof(PetscInt),&ialloc);CHKERRQ(ierr);
574    ihlp1 = ialloc;
575    ihlp2 = icol;
576 
577    if (val) {
578       ierr = PetscMalloc(nnz*sizeof(PetscScalar),&valloc);CHKERRQ(ierr);
579       vhlp1 = valloc;
580       vhlp2 = val;
581    }
582 
583 
584    /* Sorted array chunks are first 1 long, and increase until they are the complete array */
585    for (istep=1; istep<nnz; istep*=2) {
586      /*
587        Combine sorted parts
588             istart:istart+istep-1 and istart+istep-1:istart+2*istep-1
589        of ihlp2 and vhlp2
590 
591        into one sorted part
592             istart:istart+2*istep-1
593        of ihlp1 and vhlp1
594      */
595       for (istart=0; istart<nnz; istart+=2*istep) {
596 
597          /* Set counters and bound array part endings */
598          i1=istart;        i1end = i1+istep;  if (i1end>nnz) {i1end=nnz;}
599          i2=istart+istep;  i2end = i2+istep;  if (i2end>nnz) {i2end=nnz;}
600 
601          /* Merge the two array parts */
602          if (val) {
603             for (i=istart; i<i2end; i++) {
604                if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) {
605                    ihlp1[i] = ihlp2[i1];
606                    vhlp1[i] = vhlp2[i1];
607                    i1++;
608                } else if (i2<i2end) {
609                    ihlp1[i] = ihlp2[i2];
610                    vhlp1[i] = vhlp2[i2];
611                    i2++;
612                } else {
613                    ihlp1[i] = ihlp2[i1];
614                    vhlp1[i] = vhlp2[i1];
615                    i1++;
616                }
617             }
618          } else {
619             for (i=istart; i<i2end; i++) {
620                if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) {
621                    ihlp1[i] = ihlp2[i1];
622                    i1++;
623                } else if (i2<i2end) {
624                    ihlp1[i] = ihlp2[i2];
625                    i2++;
626                } else {
627                    ihlp1[i] = ihlp2[i1];
628                    i1++;
629                }
630             }
631          }
632       }
633 
634       /* Swap the two array sets */
635       iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap;
636       vswap = vhlp2; vhlp2 = vhlp1; vhlp1 = vswap;
637    }
638 
639    /* Copy one more time in case the sorted arrays are the temporary ones */
640    if (ihlp2 != icol) {
641       for (i=0; i<nnz; i++) { icol[i] = ihlp2[i]; }
642       if (val) {
643          for (i=0; i<nnz; i++) { val[i] = vhlp2[i]; }
644       }
645    }
646 
647    ierr = PetscFree(ialloc);CHKERRQ(ierr);
648    if (val) {ierr = PetscFree(valloc);CHKERRQ(ierr);}
649    PetscFunctionReturn(0);
650 }
651 
652 /*
653   spbas_apply_reordering_rows:
654     apply the given reordering to the rows:  matrix_A = matrix_A(perm,:);
655 */
656 #undef __FUNCT__
657 #define __FUNCT__ "spbas_apply_reordering_rows"
658 PetscErrorCode spbas_apply_reordering_rows(spbas_matrix *matrix_A, const PetscInt *permutation)
659 {
660    PetscInt       i,j,ip;
661    PetscInt       nrows=matrix_A->nrows;
662    PetscInt       * row_nnz;
663    PetscInt       **icols;
664    PetscBool      do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE;
665    PetscScalar    **vals=PETSC_NULL;
666    PetscErrorCode ierr;
667 
668    PetscFunctionBegin;
669    if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern\n");
670 
671    if (do_values)  {
672      ierr = PetscMalloc(sizeof(PetscScalar*)*nrows, &vals);CHKERRQ(ierr);
673    }
674    ierr = PetscMalloc(sizeof(PetscInt)*nrows, &row_nnz);CHKERRQ(ierr);
675    ierr = PetscMalloc(sizeof(PetscInt*)*nrows, &icols);CHKERRQ(ierr);
676 
677    for (i=0; i<nrows;i++) {
678       ip = permutation[i];
679       if (do_values) {vals[i]    = matrix_A->values[ip];}
680       icols[i]   = matrix_A->icols[ip];
681       row_nnz[i] = matrix_A->row_nnz[ip];
682       for (j=0; j<row_nnz[i]; j++) { icols[i][j] += ip-i; }
683    }
684 
685    if (do_values) { ierr = PetscFree(matrix_A->values);CHKERRQ(ierr);}
686    ierr = PetscFree(matrix_A->icols);CHKERRQ(ierr);
687    ierr = PetscFree(matrix_A->row_nnz);CHKERRQ(ierr);
688 
689    if (do_values) { matrix_A->values  = vals; }
690    matrix_A->icols   = icols;
691    matrix_A->row_nnz = row_nnz;
692 
693    PetscFunctionReturn(0);
694 }
695 
696 
697 /*
698   spbas_apply_reordering_cols:
699     apply the given reordering to the columns:  matrix_A(:,perm) = matrix_A;
700 */
701 #undef __FUNCT__
702 #define __FUNCT__ "spbas_apply_reordering_cols"
703 PetscErrorCode spbas_apply_reordering_cols(spbas_matrix *matrix_A,const PetscInt *permutation)
704 {
705    PetscInt    i,j;
706    PetscInt    nrows=matrix_A->nrows;
707    PetscInt    row_nnz;
708    PetscInt    *icols;
709    PetscBool   do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE;
710    PetscScalar *vals=PETSC_NULL;
711    PetscErrorCode ierr;
712 
713    PetscFunctionBegin;
714    if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "must have diagonal offsets in pattern\n");
715 
716    for (i=0; i<nrows;i++) {
717       icols   = matrix_A->icols[i];
718       row_nnz = matrix_A->row_nnz[i];
719       if (do_values) { vals = matrix_A->values[i]; }
720 
721       for (j=0; j<row_nnz; j++)  {
722          icols[j] = permutation[i+icols[j]]-i;
723       }
724       ierr = spbas_mergesort(row_nnz, icols, vals);CHKERRQ(ierr);
725    }
726    PetscFunctionReturn(0);
727 }
728 
729 /*
730   spbas_apply_reordering:
731     apply the given reordering:  matrix_A(perm,perm) = matrix_A;
732 */
733 #undef __FUNCT__
734 #define __FUNCT__ "spbas_apply_reordering"
735 PetscErrorCode spbas_apply_reordering(spbas_matrix *matrix_A, const PetscInt *permutation, const PetscInt * inv_perm)
736 {
737    PetscErrorCode ierr;
738    PetscFunctionBegin;
739    ierr = spbas_apply_reordering_rows(matrix_A, inv_perm);CHKERRQ(ierr);
740    ierr = spbas_apply_reordering_cols(matrix_A, permutation);CHKERRQ(ierr);
741    PetscFunctionReturn(0);
742 }
743 
744 #undef __FUNCT__
745 #define __FUNCT__ "spbas_pattern_only"
746 PetscErrorCode spbas_pattern_only(PetscInt nrows, PetscInt ncols, PetscInt *ai, PetscInt *aj, spbas_matrix * result)
747 {
748    spbas_matrix   retval;
749    PetscInt       i, j, i0, r_nnz;
750    PetscErrorCode ierr;
751 
752    PetscFunctionBegin;
753    /* Copy input values */
754    retval.nrows = nrows;
755    retval.ncols = ncols;
756    retval.nnz   = ai[nrows];
757 
758    retval.block_data   = PETSC_TRUE;
759    retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS;
760 
761    /* Allocate output matrix */
762    ierr =  spbas_allocate_pattern(&retval, PETSC_FALSE);CHKERRQ(ierr);
763    for (i=0; i<nrows; i++)  {retval.row_nnz[i] = ai[i+1]-ai[i];}
764    ierr =  spbas_allocate_data(&retval);CHKERRQ(ierr);
765    /* Copy the structure */
766    for (i = 0; i<retval.nrows; i++)  {
767       i0 = ai[i];
768       r_nnz = ai[i+1]-i0;
769 
770       for (j=0; j<r_nnz; j++) {
771          retval.icols[i][j]  = aj[i0+j]-i;
772       }
773    }
774    *result = retval;
775    PetscFunctionReturn(0);
776 }
777 
778 
779 /*
780    spbas_mark_row_power:
781       Mark the columns in row 'row' which are nonzero in
782           matrix^2log(marker).
783 */
784 #undef __FUNCT__
785 #define __FUNCT__ "spbas_mark_row_power"
786 PetscErrorCode spbas_mark_row_power(
787                                     PetscInt *iwork,             /* marker-vector */
788                                     PetscInt row,                /* row for which the columns are marked */
789                                     spbas_matrix * in_matrix, /* matrix for which the power is being  calculated */
790                                     PetscInt marker,             /* marker-value: 2^power */
791                                     PetscInt minmrk,            /* lower bound for marked points */
792                                     PetscInt maxmrk)            /* upper bound for marked points */
793 {
794    PetscErrorCode ierr;
795    PetscInt       i,j, nnz;
796 
797    PetscFunctionBegin;
798    nnz = in_matrix->row_nnz[row];
799 
800    /* For higher powers, call this function recursively */
801    if (marker>1) {
802       for (i=0; i<nnz;i++) {
803          j = row + in_matrix->icols[row][i];
804          if (minmrk<=j && j<maxmrk && iwork[j] < marker) {
805             ierr =  spbas_mark_row_power(iwork, row + in_matrix->icols[row][i],in_matrix, marker/2,minmrk,maxmrk);CHKERRQ(ierr);
806             iwork[j] |= marker;
807          }
808       }
809    } else  {
810      /*  Mark the columns reached */
811       for (i=0; i<nnz;i++)  {
812          j = row + in_matrix->icols[row][i];
813          if (minmrk<=j && j<maxmrk) { iwork[j] |= 1; }
814       }
815    }
816    PetscFunctionReturn(0);
817 }
818 
819 
820 /*
821    spbas_power
822       Calculate sparseness patterns for incomplete Cholesky decompositions
823       of a given order: (almost) all nonzeros of the matrix^(order+1) which
824       are inside the band width are found and stored in the output sparseness
825       pattern.
826 */
827 #undef __FUNCT__
828 #define __FUNCT__ "spbas_power"
829 PetscErrorCode spbas_power (spbas_matrix in_matrix,PetscInt power, spbas_matrix * result)
830 {
831    spbas_matrix   retval;
832    PetscInt       nrows = in_matrix.nrows;
833    PetscInt       ncols = in_matrix.ncols;
834    PetscInt       i, j, kend;
835    PetscInt       nnz, inz;
836    PetscInt       *iwork;
837    PetscInt       marker;
838    PetscInt       maxmrk=0;
839    PetscErrorCode ierr;
840 
841    PetscFunctionBegin;
842    if (in_matrix.col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern\n");
843    if (ncols != nrows) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Dimension error\n");
844    if (in_matrix.values) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input array must be sparseness pattern (no values)");
845    if (power<=0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Power must be 1 or up");
846 
847    /* Copy input values*/
848    retval.nrows = ncols;
849    retval.ncols = nrows;
850    retval.nnz   = 0;
851    retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS;
852    retval.block_data = PETSC_FALSE;
853 
854    /* Allocate sparseness pattern */
855    ierr =  spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
856 
857    /* Allocate marker array */
858    ierr = PetscMalloc(nrows * sizeof(PetscInt), &iwork);CHKERRQ(ierr);
859 
860    /* Erase the pattern for this row */
861    ierr = PetscMemzero((void *) iwork, retval.nrows*sizeof(PetscInt));CHKERRQ(ierr);
862 
863    /* Calculate marker values */
864    marker = 1; for (i=1; i<power; i++) {marker*=2;}
865 
866    for (i=0; i<nrows; i++)  {
867      /* Calculate the pattern for each row */
868 
869       nnz    = in_matrix.row_nnz[i];
870       kend   = i+in_matrix.icols[i][nnz-1];
871       if (maxmrk<=kend) {maxmrk=kend+1;}
872       ierr =  spbas_mark_row_power(iwork, i, &in_matrix, marker, i, maxmrk);CHKERRQ(ierr);
873 
874       /* Count the columns*/
875       nnz = 0;
876       for (j=i; j<maxmrk; j++) {nnz+= (iwork[j]!=0);}
877 
878       /* Allocate the column indices */
879       retval.row_nnz[i] = nnz;
880       ierr = PetscMalloc(nnz*sizeof(PetscInt),&retval.icols[i]);CHKERRQ(ierr);
881 
882       /* Administrate the column indices */
883       inz = 0;
884       for (j=i; j<maxmrk; j++) {
885          if (iwork[j])  {
886             retval.icols[i][inz] = j-i;
887             inz++;
888             iwork[j]=0;
889          }
890       }
891       retval.nnz += nnz;
892    };
893    ierr = PetscFree(iwork);CHKERRQ(ierr);
894    *result = retval;
895    PetscFunctionReturn(0);
896 }
897 
898 
899 
900 /*
901    spbas_keep_upper:
902       remove the lower part of the matrix: keep the upper part
903 */
904 #undef __FUNCT__
905 #define __FUNCT__ "spbas_keep_upper"
906 PetscErrorCode spbas_keep_upper(spbas_matrix * inout_matrix)
907 {
908    PetscInt i, j;
909    PetscInt jstart;
910 
911    PetscFunctionBegin;
912    if (inout_matrix->block_data) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Not yet for block data matrices\n");
913    for (i=0; i<inout_matrix->nrows; i++)  {
914        for (jstart=0; (jstart<inout_matrix->row_nnz[i]) && (inout_matrix->icols[i][jstart]<0); jstart++) {}
915        if (jstart>0) {
916           for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) {
917              inout_matrix->icols[i][j] =  inout_matrix->icols[i][j+jstart];
918           }
919 
920           if (inout_matrix->values) {
921              for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) {
922                 inout_matrix->values[i][j] = inout_matrix->values[i][j+jstart];
923              }
924           }
925 
926           inout_matrix->row_nnz[i] -= jstart;
927 
928           inout_matrix->icols[i] = (PetscInt*) realloc((void*) inout_matrix->icols[i], inout_matrix->row_nnz[i]*sizeof(PetscInt));
929 
930           if (inout_matrix->values) {
931              inout_matrix->values[i] = (PetscScalar*) realloc((void*) inout_matrix->values[i], inout_matrix->row_nnz[i]*sizeof(PetscScalar));
932           }
933           inout_matrix->nnz -= jstart;
934        }
935    }
936    PetscFunctionReturn(0);
937 }
938 
939