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