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