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