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