xref: /petsc/src/mat/impls/aij/seq/bas/spbas.c (revision 580bdb303e1ee3b1222b2042810b4c26340259c6)
1c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas.h>
3b2f1ab58SBarry Smith 
4845006b9SBarry Smith /*MC
52692d6eeSBarry Smith   MATSOLVERBAS -  Provides ICC(k) with drop tolerance
6845006b9SBarry Smith 
7845006b9SBarry Smith   Works with MATAIJ  matrices
8845006b9SBarry Smith 
9845006b9SBarry Smith   Options Database Keys:
100053dfc8SBarry Smith + -pc_factor_levels <l> - number of levels of fill
110053dfc8SBarry Smith - -pc_factor_drop_tolerance - is not currently hooked up to do anything
12845006b9SBarry Smith 
130053dfc8SBarry Smith   Level: intermediate
14845006b9SBarry Smith 
15845006b9SBarry Smith    Contributed by: Bas van 't Hof
16845006b9SBarry Smith 
1795452b02SPatrick Sanan    Notes:
1895452b02SPatrick Sanan     Since this currently hooked up to use drop tolerance it should produce the same factors and hence convergence as the PETSc ICC, for higher
190053dfc8SBarry Smith      levels of fill it does not. This needs to be investigated. Unless you are interested in drop tolerance ICC and willing to work through the code
200053dfc8SBarry Smith      we recommend not using this functionality.
210053dfc8SBarry Smith 
223ca39a21SBarry Smith .seealso: PCFactorSetMatSolverType(), MatSolverType, PCFactorSetLevels(), PCFactorSetDropTolerance()
23845006b9SBarry Smith 
24845006b9SBarry Smith M*/
259ccc4a9bSBarry Smith 
26b2f1ab58SBarry Smith /*
27b2f1ab58SBarry Smith   spbas_memory_requirement:
28b2f1ab58SBarry Smith     Calculate the number of bytes needed to store tha matrix
29b2f1ab58SBarry Smith */
303dfa2556SBarry Smith size_t spbas_memory_requirement(spbas_matrix matrix)
31b2f1ab58SBarry Smith {
323dfa2556SBarry Smith   size_t memreq = 6 * sizeof(PetscInt)  + /* nrows, ncols, nnz, n_alloc_icol, n_alloc_val, col_idx_type */
33ace3abfcSBarry Smith                     sizeof(PetscBool)               + /* block_data */
34c328eeadSBarry Smith                     sizeof(PetscScalar**)           + /* values */
35c328eeadSBarry Smith                     sizeof(PetscScalar*)            + /* alloc_val */
363dfa2556SBarry Smith                     2 * sizeof(PetscInt**)          + /* icols, icols0 */
37c328eeadSBarry Smith                     2 * sizeof(PetscInt*)           + /* row_nnz, alloc_icol */
38c328eeadSBarry Smith                     matrix.nrows * sizeof(PetscInt) + /* row_nnz[*] */
39c328eeadSBarry Smith                     matrix.nrows * sizeof(PetscInt*); /* icols[*] */
40b2f1ab58SBarry Smith 
41c328eeadSBarry Smith   /* icol0[*] */
422205254eSKarl Rupp   if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) memreq += matrix.nrows * sizeof(PetscInt);
43b2f1ab58SBarry Smith 
44c328eeadSBarry Smith   /* icols[*][*] */
452205254eSKarl Rupp   if (matrix.block_data) memreq += matrix.n_alloc_icol * sizeof(PetscInt);
462205254eSKarl Rupp   else memreq += matrix.nnz * sizeof(PetscInt);
47b2f1ab58SBarry Smith 
484e1ff37aSBarry Smith   if (matrix.values) {
49c328eeadSBarry Smith     memreq += matrix.nrows * sizeof(PetscScalar*); /* values[*] */
50c328eeadSBarry Smith     /* values[*][*] */
512205254eSKarl Rupp     if (matrix.block_data) memreq += matrix.n_alloc_val  * sizeof(PetscScalar);
522205254eSKarl Rupp     else memreq += matrix.nnz * sizeof(PetscScalar);
53b2f1ab58SBarry Smith   }
54b2f1ab58SBarry Smith   return memreq;
55b2f1ab58SBarry Smith }
56b2f1ab58SBarry Smith 
57b2f1ab58SBarry Smith /*
58b2f1ab58SBarry Smith   spbas_allocate_pattern:
59b2f1ab58SBarry Smith     allocate the pattern arrays row_nnz, icols and optionally values
60b2f1ab58SBarry Smith */
61ace3abfcSBarry Smith PetscErrorCode spbas_allocate_pattern(spbas_matrix * result, PetscBool do_values)
62b2f1ab58SBarry Smith {
63b2f1ab58SBarry Smith   PetscErrorCode ierr;
64b2f1ab58SBarry Smith   PetscInt       nrows        = result->nrows;
65b2f1ab58SBarry Smith   PetscInt       col_idx_type = result->col_idx_type;
66b2f1ab58SBarry Smith 
67b2f1ab58SBarry Smith   PetscFunctionBegin;
68c328eeadSBarry Smith   /* Allocate sparseness pattern */
69785e854fSJed Brown   ierr = PetscMalloc1(nrows,&result->row_nnz);CHKERRQ(ierr);
70785e854fSJed Brown   ierr = PetscMalloc1(nrows,&result->icols);CHKERRQ(ierr);
71b2f1ab58SBarry Smith 
72c328eeadSBarry Smith   /* If offsets are given wrt an array, create array */
734e1ff37aSBarry Smith   if (col_idx_type == SPBAS_OFFSET_ARRAY) {
74785e854fSJed Brown     ierr = PetscMalloc1(nrows,&result->icol0);CHKERRQ(ierr);
754e1ff37aSBarry Smith   } else  {
760298fd71SBarry Smith     result->icol0 = NULL;
77b2f1ab58SBarry Smith   }
78b2f1ab58SBarry Smith 
79c328eeadSBarry Smith   /* If values are given, allocate values array */
804e1ff37aSBarry Smith   if (do_values)  {
81785e854fSJed Brown     ierr = PetscMalloc1(nrows,&result->values);CHKERRQ(ierr);
824e1ff37aSBarry Smith   } else {
830298fd71SBarry Smith     result->values = NULL;
84b2f1ab58SBarry Smith   }
85b2f1ab58SBarry Smith   PetscFunctionReturn(0);
86b2f1ab58SBarry Smith }
87b2f1ab58SBarry Smith 
88b2f1ab58SBarry Smith /*
89b2f1ab58SBarry Smith spbas_allocate_data:
90b2f1ab58SBarry Smith    in case of block_data:
91b2f1ab58SBarry Smith        Allocate the data arrays alloc_icol and optionally alloc_val,
92b2f1ab58SBarry Smith        set appropriate pointers from icols and values;
93b2f1ab58SBarry Smith    in case of !block_data:
94b2f1ab58SBarry Smith        Allocate the arrays icols[i] and optionally values[i]
95b2f1ab58SBarry Smith */
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;
1036c4ed002SBarry Smith   PetscBool      do_values  = (result->values) ? PETSC_TRUE : PETSC_FALSE;
104ace3abfcSBarry Smith   PetscBool      block_data = result->block_data;
105b2f1ab58SBarry Smith 
106b2f1ab58SBarry Smith   PetscFunctionBegin;
1074e1ff37aSBarry Smith   if (block_data) {
108c328eeadSBarry Smith     /* Allocate the column number array and point to it */
109b2f1ab58SBarry Smith     result->n_alloc_icol = nnz;
1102205254eSKarl Rupp 
111785e854fSJed Brown     ierr = PetscMalloc1(nnz, &result->alloc_icol);CHKERRQ(ierr);
1122205254eSKarl Rupp 
113b2f1ab58SBarry Smith     result->icols[0] = result->alloc_icol;
1144e1ff37aSBarry Smith     for (i=1; i<nrows; i++)  {
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 */
1194e1ff37aSBarry Smith     if (do_values) {
120b2f1ab58SBarry Smith       result->n_alloc_val = nnz;
1212205254eSKarl Rupp 
122785e854fSJed Brown       ierr = PetscMalloc1(nnz, &result->alloc_val);CHKERRQ(ierr);
1232205254eSKarl Rupp 
124b2f1ab58SBarry Smith       result->values[0] = result->alloc_val;
1254e1ff37aSBarry Smith       for (i=1; i<nrows; i++) {
126b2f1ab58SBarry Smith         result->values[i] = result->values[i-1] + result->row_nnz[i-1];
127b2f1ab58SBarry Smith       }
128b2f1ab58SBarry Smith     }
1294e1ff37aSBarry Smith   } else {
1304e1ff37aSBarry Smith     for (i=0; i<nrows; i++)   {
131b2f1ab58SBarry Smith       r_nnz = result->row_nnz[i];
132785e854fSJed Brown       ierr  = PetscMalloc1(r_nnz, &result->icols[i]);CHKERRQ(ierr);
133b2f1ab58SBarry Smith     }
1344e1ff37aSBarry Smith     if (do_values) {
1354e1ff37aSBarry Smith       for (i=0; i<nrows; i++)  {
136b2f1ab58SBarry Smith         r_nnz = result->row_nnz[i];
137785e854fSJed Brown         ierr  = PetscMalloc1(r_nnz, &result->values[i]);CHKERRQ(ierr);
138b2f1ab58SBarry Smith       }
139b2f1ab58SBarry Smith     }
140b2f1ab58SBarry Smith   }
141b2f1ab58SBarry Smith   PetscFunctionReturn(0);
142b2f1ab58SBarry Smith }
143b2f1ab58SBarry Smith 
144b2f1ab58SBarry Smith /*
145b2f1ab58SBarry Smith   spbas_row_order_icol
146b2f1ab58SBarry Smith      determine if row i1 should come
147b2f1ab58SBarry Smith        + before row i2 in the sorted rows (return -1),
148b2f1ab58SBarry Smith        + after (return 1)
149b2f1ab58SBarry Smith        + is identical (return 0).
150b2f1ab58SBarry Smith */
151b2f1ab58SBarry Smith int spbas_row_order_icol(PetscInt i1, PetscInt i2, PetscInt *irow_in, PetscInt *icol_in,PetscInt col_idx_type)
152b2f1ab58SBarry Smith {
153b2f1ab58SBarry Smith   PetscInt j;
154b2f1ab58SBarry Smith   PetscInt nnz1    = irow_in[i1+1] - irow_in[i1];
155b2f1ab58SBarry Smith   PetscInt nnz2    = irow_in[i2+1] - irow_in[i2];
156b2f1ab58SBarry Smith   PetscInt * icol1 = &icol_in[irow_in[i1]];
157b2f1ab58SBarry Smith   PetscInt * icol2 = &icol_in[irow_in[i2]];
158b2f1ab58SBarry Smith 
1592205254eSKarl Rupp   if (nnz1<nnz2) return -1;
1602205254eSKarl Rupp   if (nnz1>nnz2) return 1;
161b2f1ab58SBarry Smith 
1624e1ff37aSBarry Smith   if (col_idx_type == SPBAS_COLUMN_NUMBERS) {
1634e1ff37aSBarry Smith     for (j=0; j<nnz1; j++) {
1642205254eSKarl Rupp       if (icol1[j]< icol2[j]) return -1;
1652205254eSKarl Rupp       if (icol1[j]> icol2[j]) return 1;
166b2f1ab58SBarry Smith     }
1674e1ff37aSBarry Smith   } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {
1684e1ff37aSBarry Smith     for (j=0; j<nnz1; j++) {
1692205254eSKarl Rupp       if (icol1[j]-i1< icol2[j]-i2) return -1;
1702205254eSKarl Rupp       if (icol1[j]-i1> icol2[j]-i2) return 1;
171b2f1ab58SBarry Smith     }
1724e1ff37aSBarry Smith   } else if (col_idx_type == SPBAS_OFFSET_ARRAY) {
1734e1ff37aSBarry Smith     for (j=1; j<nnz1; j++) {
1742205254eSKarl Rupp       if (icol1[j]-icol1[0] < icol2[j]-icol2[0]) return -1;
1752205254eSKarl Rupp       if (icol1[j]-icol1[0] > icol2[j]-icol2[0]) return 1;
176b2f1ab58SBarry Smith     }
177b2f1ab58SBarry Smith   }
178b2f1ab58SBarry Smith   return 0;
179b2f1ab58SBarry Smith }
180b2f1ab58SBarry Smith 
181b2f1ab58SBarry Smith /*
182b2f1ab58SBarry Smith   spbas_mergesort_icols:
183b2f1ab58SBarry Smith     return a sorting of the rows in which identical sparseness patterns are
184b2f1ab58SBarry Smith     next to each other
185b2f1ab58SBarry Smith */
186b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort_icols(PetscInt nrows, PetscInt * irow_in, PetscInt * icol_in,PetscInt col_idx_type, PetscInt *isort)
187b2f1ab58SBarry Smith {
188b2f1ab58SBarry Smith   PetscErrorCode ierr;
189c328eeadSBarry Smith   PetscInt       istep;       /* Chunk-sizes of already sorted parts of arrays */
190c328eeadSBarry Smith   PetscInt       i, i1, i2;   /* Loop counters for (partly) sorted arrays */
191c328eeadSBarry Smith   PetscInt       istart, i1end, i2end; /* start of newly sorted array part, end of both  parts */
192c328eeadSBarry Smith   PetscInt       *ialloc;     /* Allocated arrays */
193c328eeadSBarry Smith   PetscInt       *iswap;      /* auxiliary pointers for swapping */
194c328eeadSBarry Smith   PetscInt       *ihlp1;      /* Pointers to new version of arrays, */
195c328eeadSBarry Smith   PetscInt       *ihlp2;      /* Pointers to previous version of arrays, */
196b2f1ab58SBarry Smith 
197b2f1ab58SBarry Smith   PetscFunctionBegin;
198785e854fSJed Brown   ierr = PetscMalloc1(nrows,&ialloc);CHKERRQ(ierr);
199b2f1ab58SBarry Smith 
200b2f1ab58SBarry Smith   ihlp1 = ialloc;
201b2f1ab58SBarry Smith   ihlp2 = isort;
202b2f1ab58SBarry Smith 
203c328eeadSBarry Smith   /* Sorted array chunks are first 1 long, and increase until they are the complete array */
2044e1ff37aSBarry Smith   for (istep=1; istep<nrows; istep*=2) {
205c328eeadSBarry Smith     /*
206c328eeadSBarry Smith       Combine sorted parts
207c328eeadSBarry Smith           istart:istart+istep-1 and istart+istep-1:istart+2*istep-1
208c328eeadSBarry Smith       of ihlp2 and vhlp2
209c328eeadSBarry Smith 
210c328eeadSBarry Smith       into one sorted part
211c328eeadSBarry Smith           istart:istart+2*istep-1
212c328eeadSBarry Smith       of ihlp1 and vhlp1
213c328eeadSBarry Smith     */
2144e1ff37aSBarry Smith     for (istart=0; istart<nrows; istart+=2*istep) {
215c328eeadSBarry Smith       /* Set counters and bound array part endings */
2162205254eSKarl Rupp       i1=istart;        i1end = i1+istep;  if (i1end>nrows) i1end=nrows;
2172205254eSKarl Rupp       i2=istart+istep;  i2end = i2+istep;  if (i2end>nrows) i2end=nrows;
218b2f1ab58SBarry Smith 
219c328eeadSBarry Smith       /* Merge the two array parts */
2204e1ff37aSBarry Smith       for (i=istart; i<i2end; i++)  {
2214e1ff37aSBarry Smith         if (i1<i1end && i2<i2end && spbas_row_order_icol(ihlp2[i1], ihlp2[i2], irow_in, icol_in, col_idx_type) < 0) {
222b2f1ab58SBarry Smith           ihlp1[i] = ihlp2[i1];
223b2f1ab58SBarry Smith           i1++;
2244e1ff37aSBarry Smith         } else if (i2<i2end) {
225b2f1ab58SBarry Smith           ihlp1[i] = ihlp2[i2];
226b2f1ab58SBarry Smith           i2++;
2274e1ff37aSBarry Smith         } else {
228b2f1ab58SBarry Smith           ihlp1[i] = ihlp2[i1];
229b2f1ab58SBarry Smith           i1++;
230b2f1ab58SBarry Smith         }
231b2f1ab58SBarry Smith       }
232b2f1ab58SBarry Smith     }
233b2f1ab58SBarry Smith 
234c328eeadSBarry Smith     /* Swap the two array sets */
235b2f1ab58SBarry Smith     iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap;
236b2f1ab58SBarry Smith   }
237b2f1ab58SBarry Smith 
238c328eeadSBarry Smith   /* Copy one more time in case the sorted arrays are the temporary ones */
2394e1ff37aSBarry Smith   if (ihlp2 != isort) {
2402205254eSKarl Rupp     for (i=0; i<nrows; i++) isort[i] = ihlp2[i];
241b2f1ab58SBarry Smith   }
242b2f1ab58SBarry Smith   ierr = PetscFree(ialloc);CHKERRQ(ierr);
243b2f1ab58SBarry Smith   PetscFunctionReturn(0);
244b2f1ab58SBarry Smith }
245b2f1ab58SBarry Smith 
246b2f1ab58SBarry Smith 
247b2f1ab58SBarry Smith 
248b2f1ab58SBarry Smith /*
249b2f1ab58SBarry Smith   spbas_compress_pattern:
250b2f1ab58SBarry Smith      calculate a compressed sparseness pattern for a sparseness pattern
251b2f1ab58SBarry Smith      given in compressed row storage. The compressed sparseness pattern may
252b2f1ab58SBarry Smith      require (much) less memory.
253b2f1ab58SBarry Smith */
2544e1ff37aSBarry Smith PetscErrorCode spbas_compress_pattern(PetscInt *irow_in, PetscInt *icol_in, PetscInt nrows, PetscInt ncols, PetscInt col_idx_type, spbas_matrix *B,PetscReal *mem_reduction)
255b2f1ab58SBarry Smith {
256b2f1ab58SBarry Smith   PetscInt        nnz      = irow_in[nrows];
2573dfa2556SBarry Smith   size_t          mem_orig = (nrows + nnz) * sizeof(PetscInt);
2583dfa2556SBarry Smith   size_t          mem_compressed;
259b2f1ab58SBarry Smith   PetscErrorCode  ierr;
260b2f1ab58SBarry Smith   PetscInt        *isort;
261b2f1ab58SBarry Smith   PetscInt        *icols;
262b2f1ab58SBarry Smith   PetscInt        row_nnz;
263b2f1ab58SBarry Smith   PetscInt        *ipoint;
264ace3abfcSBarry Smith   PetscBool       *used;
265b2f1ab58SBarry Smith   PetscInt        ptr;
266b2f1ab58SBarry Smith   PetscInt        i,j;
267ace3abfcSBarry Smith   const PetscBool no_values = PETSC_FALSE;
268b2f1ab58SBarry Smith 
269b2f1ab58SBarry Smith   PetscFunctionBegin;
270c328eeadSBarry Smith   /* Allocate the structure of the new matrix */
271b2f1ab58SBarry Smith   B->nrows        = nrows;
272b2f1ab58SBarry Smith   B->ncols        = ncols;
273b2f1ab58SBarry Smith   B->nnz          = nnz;
274b2f1ab58SBarry Smith   B->col_idx_type = col_idx_type;
275b2f1ab58SBarry Smith   B->block_data   = PETSC_TRUE;
2762205254eSKarl Rupp 
277b2f1ab58SBarry Smith   ierr = spbas_allocate_pattern(B, no_values);CHKERRQ(ierr);
278b2f1ab58SBarry Smith 
279c328eeadSBarry Smith   /* When using an offset array, set it */
2804e1ff37aSBarry Smith   if (col_idx_type==SPBAS_OFFSET_ARRAY)  {
2812205254eSKarl Rupp     for (i=0; i<nrows; i++) B->icol0[i] = icol_in[irow_in[i]];
282b2f1ab58SBarry Smith   }
283b2f1ab58SBarry Smith 
284c328eeadSBarry Smith   /* Allocate the ordering for the rows */
285785e854fSJed Brown   ierr = PetscMalloc1(nrows,&isort);CHKERRQ(ierr);
286785e854fSJed Brown   ierr = PetscMalloc1(nrows,&ipoint);CHKERRQ(ierr);
287*580bdb30SBarry Smith   ierr = PetscCalloc1(nrows,&used);CHKERRQ(ierr);
288b2f1ab58SBarry Smith 
2894e1ff37aSBarry Smith   for (i = 0; i<nrows; i++)  {
290b2f1ab58SBarry Smith     B->row_nnz[i] = irow_in[i+1]-irow_in[i];
291b2f1ab58SBarry Smith     isort[i]      = i;
292b2f1ab58SBarry Smith     ipoint[i]     = i;
293b2f1ab58SBarry Smith   }
294b2f1ab58SBarry Smith 
295c328eeadSBarry Smith   /* Sort the rows so that identical columns will be next to each other */
296b2f1ab58SBarry Smith   ierr = spbas_mergesort_icols(nrows, irow_in, icol_in, col_idx_type, isort);CHKERRQ(ierr);
2970298fd71SBarry Smith   ierr = PetscInfo(NULL,"Rows have been sorted for patterns\n");CHKERRQ(ierr);
298b2f1ab58SBarry Smith 
299c328eeadSBarry Smith   /* Replace identical rows with the first one in the list */
3004e1ff37aSBarry Smith   for (i=1; i<nrows; i++) {
3014e1ff37aSBarry Smith     if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) {
302b2f1ab58SBarry Smith       ipoint[isort[i]] = ipoint[isort[i-1]];
303b2f1ab58SBarry Smith     }
304b2f1ab58SBarry Smith   }
305b2f1ab58SBarry Smith 
306c328eeadSBarry Smith   /* Collect the rows which are used*/
3072205254eSKarl Rupp   for (i=0; i<nrows; i++) used[ipoint[i]] = PETSC_TRUE;
308b2f1ab58SBarry Smith 
309c328eeadSBarry Smith   /* Calculate needed memory */
310b2f1ab58SBarry Smith   B->n_alloc_icol = 0;
3114e1ff37aSBarry Smith   for (i=0; i<nrows; i++)  {
3122205254eSKarl Rupp     if (used[i]) B->n_alloc_icol += B->row_nnz[i];
313b2f1ab58SBarry Smith   }
314785e854fSJed Brown   ierr = PetscMalloc1(B->n_alloc_icol,&B->alloc_icol);CHKERRQ(ierr);
315b2f1ab58SBarry Smith 
316c328eeadSBarry Smith   /* Fill in the diagonal offsets for the rows which store their own data */
317b2f1ab58SBarry Smith   ptr = 0;
3184e1ff37aSBarry Smith   for (i=0; i<B->nrows; i++) {
3194e1ff37aSBarry Smith     if (used[i]) {
320b2f1ab58SBarry Smith       B->icols[i] = &B->alloc_icol[ptr];
321b2f1ab58SBarry Smith       icols = &icol_in[irow_in[i]];
322b2f1ab58SBarry Smith       row_nnz = B->row_nnz[i];
3234e1ff37aSBarry Smith       if (col_idx_type == SPBAS_COLUMN_NUMBERS) {
3244e1ff37aSBarry Smith         for (j=0; j<row_nnz; j++) {
325b2f1ab58SBarry Smith           B->icols[i][j] = icols[j];
326b2f1ab58SBarry Smith         }
3274e1ff37aSBarry Smith       } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {
3284e1ff37aSBarry Smith         for (j=0; j<row_nnz; j++) {
329b2f1ab58SBarry Smith           B->icols[i][j] = icols[j]-i;
330b2f1ab58SBarry Smith         }
3314e1ff37aSBarry Smith       } else if (col_idx_type == SPBAS_OFFSET_ARRAY) {
3324e1ff37aSBarry Smith         for (j=0; j<row_nnz; j++) {
333b2f1ab58SBarry Smith           B->icols[i][j] = icols[j]-icols[0];
334b2f1ab58SBarry Smith         }
335b2f1ab58SBarry Smith       }
336b2f1ab58SBarry Smith       ptr += B->row_nnz[i];
337b2f1ab58SBarry Smith     }
338b2f1ab58SBarry Smith   }
339b2f1ab58SBarry Smith 
340c328eeadSBarry Smith   /* Point to the right places for all data */
3414e1ff37aSBarry Smith   for (i=0; i<nrows; i++) {
342b2f1ab58SBarry Smith     B->icols[i] = B->icols[ipoint[i]];
343b2f1ab58SBarry Smith   }
3440298fd71SBarry Smith   ierr = PetscInfo(NULL,"Row patterns have been compressed\n");CHKERRQ(ierr);
34557622a8eSBarry Smith   ierr = PetscInfo1(NULL,"         (%g nonzeros per row)\n", (double) ((PetscReal) nnz / (PetscReal) nrows));CHKERRQ(ierr);
346b2f1ab58SBarry Smith 
347b2f1ab58SBarry Smith   ierr=PetscFree(isort);CHKERRQ(ierr);
348b2f1ab58SBarry Smith   ierr=PetscFree(used);CHKERRQ(ierr);
349b2f1ab58SBarry Smith   ierr=PetscFree(ipoint);CHKERRQ(ierr);
350b2f1ab58SBarry Smith 
351b2f1ab58SBarry Smith   mem_compressed = spbas_memory_requirement(*B);
3524e1ff37aSBarry Smith   *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig;
353b2f1ab58SBarry Smith   PetscFunctionReturn(0);
354b2f1ab58SBarry Smith }
355b2f1ab58SBarry Smith 
356b2f1ab58SBarry Smith /*
357b2f1ab58SBarry Smith    spbas_incomplete_cholesky
358b2f1ab58SBarry Smith        Incomplete Cholesky decomposition
359b2f1ab58SBarry Smith */
360c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas_cholesky.h>
361b2f1ab58SBarry Smith 
362b2f1ab58SBarry Smith /*
363b2f1ab58SBarry Smith   spbas_delete : de-allocate the arrays owned by this matrix
364b2f1ab58SBarry Smith */
365b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix matrix)
366b2f1ab58SBarry Smith {
367b2f1ab58SBarry Smith   PetscInt       i;
368b2f1ab58SBarry Smith   PetscErrorCode ierr;
3699ccc4a9bSBarry Smith 
370b2f1ab58SBarry Smith   PetscFunctionBegin;
3719ccc4a9bSBarry Smith   if (matrix.block_data) {
372cb9801acSJed Brown     ierr=PetscFree(matrix.alloc_icol);CHKERRQ(ierr);
373b2f1ab58SBarry Smith     if (matrix.values) {ierr=PetscFree(matrix.alloc_val);CHKERRQ(ierr);}
3749ccc4a9bSBarry Smith   } else {
3759ccc4a9bSBarry Smith     for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.icols[i]);CHKERRQ(ierr);}
376b2f1ab58SBarry Smith     ierr = PetscFree(matrix.icols);CHKERRQ(ierr);
3779ccc4a9bSBarry Smith     if (matrix.values) {
3789ccc4a9bSBarry Smith       for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.values[i]);CHKERRQ(ierr);}
379b2f1ab58SBarry Smith     }
380b2f1ab58SBarry Smith   }
381b2f1ab58SBarry Smith 
382b2f1ab58SBarry Smith   ierr=PetscFree(matrix.row_nnz);CHKERRQ(ierr);
383b2f1ab58SBarry Smith   ierr=PetscFree(matrix.icols);CHKERRQ(ierr);
3849ccc4a9bSBarry Smith   if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) {ierr=PetscFree(matrix.icol0);CHKERRQ(ierr);}
385c31cb41cSBarry Smith   ierr=PetscFree(matrix.values);CHKERRQ(ierr);
386b2f1ab58SBarry Smith   PetscFunctionReturn(0);
387b2f1ab58SBarry Smith }
388b2f1ab58SBarry Smith 
389b2f1ab58SBarry Smith /*
390b2f1ab58SBarry Smith spbas_matrix_to_crs:
391b2f1ab58SBarry Smith    Convert an spbas_matrix to compessed row storage
392b2f1ab58SBarry Smith */
393b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out)
394b2f1ab58SBarry Smith {
395b2f1ab58SBarry Smith   PetscInt       nrows = matrix_A.nrows;
396b2f1ab58SBarry Smith   PetscInt       nnz   = matrix_A.nnz;
397b2f1ab58SBarry Smith   PetscInt       i,j,r_nnz,i0;
398b2f1ab58SBarry Smith   PetscInt       *irow;
399b2f1ab58SBarry Smith   PetscInt       *icol;
400b2f1ab58SBarry Smith   PetscInt       *icol_A;
401b2f1ab58SBarry Smith   MatScalar      *val;
402b2f1ab58SBarry Smith   PetscScalar    *val_A;
403b2f1ab58SBarry Smith   PetscInt       col_idx_type = matrix_A.col_idx_type;
404ace3abfcSBarry Smith   PetscBool      do_values    = matrix_A.values ? PETSC_TRUE : PETSC_FALSE;
405b2f1ab58SBarry Smith   PetscErrorCode ierr;
406b2f1ab58SBarry Smith 
407b2f1ab58SBarry Smith   PetscFunctionBegin;
408854ce69bSBarry Smith   ierr      = PetscMalloc1(nrows+1, &irow);CHKERRQ(ierr);
409854ce69bSBarry Smith   ierr      = PetscMalloc1(nnz, &icol);CHKERRQ(ierr);
4109ccc4a9bSBarry Smith   *icol_out = icol;
4119ccc4a9bSBarry Smith   *irow_out = irow;
4129ccc4a9bSBarry Smith   if (do_values) {
413854ce69bSBarry Smith     ierr     = PetscMalloc1(nnz, &val);CHKERRQ(ierr);
414b2f1ab58SBarry Smith     *val_out = val; *icol_out = icol; *irow_out=irow;
415b2f1ab58SBarry Smith   }
416b2f1ab58SBarry Smith 
417b2f1ab58SBarry Smith   irow[0]=0;
4189ccc4a9bSBarry Smith   for (i=0; i<nrows; i++) {
419b2f1ab58SBarry Smith     r_nnz     = matrix_A.row_nnz[i];
420b2f1ab58SBarry Smith     i0        = irow[i];
421b2f1ab58SBarry Smith     irow[i+1] = i0 + r_nnz;
422b2f1ab58SBarry Smith     icol_A    = matrix_A.icols[i];
423b2f1ab58SBarry Smith 
4249ccc4a9bSBarry Smith     if (do_values) {
425b2f1ab58SBarry Smith       val_A = matrix_A.values[i];
4269ccc4a9bSBarry Smith       for (j=0; j<r_nnz; j++) {
427b2f1ab58SBarry Smith         icol[i0+j] = icol_A[j];
428b2f1ab58SBarry Smith         val[i0+j]  = val_A[j];
429b2f1ab58SBarry Smith       }
4309ccc4a9bSBarry Smith     } else {
4312205254eSKarl Rupp       for (j=0; j<r_nnz; j++) icol[i0+j] = icol_A[j];
432b2f1ab58SBarry Smith     }
433b2f1ab58SBarry Smith 
4349ccc4a9bSBarry Smith     if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {
4352205254eSKarl Rupp       for (j=0; j<r_nnz; j++) icol[i0+j] += i;
4362205254eSKarl Rupp     } else if (col_idx_type == SPBAS_OFFSET_ARRAY) {
437b2f1ab58SBarry Smith       i0 = matrix_A.icol0[i];
4382205254eSKarl Rupp       for (j=0; j<r_nnz; j++) icol[i0+j] += i0;
439b2f1ab58SBarry Smith     }
440b2f1ab58SBarry Smith   }
441b2f1ab58SBarry Smith   PetscFunctionReturn(0);
442b2f1ab58SBarry Smith }
443b2f1ab58SBarry Smith 
444b2f1ab58SBarry Smith 
445b2f1ab58SBarry Smith /*
446b2f1ab58SBarry Smith     spbas_transpose
447b2f1ab58SBarry Smith        return the transpose of a matrix
448b2f1ab58SBarry Smith */
449b2f1ab58SBarry Smith PetscErrorCode spbas_transpose(spbas_matrix in_matrix, spbas_matrix * result)
450b2f1ab58SBarry Smith {
451b2f1ab58SBarry Smith   PetscInt       col_idx_type = in_matrix.col_idx_type;
452b2f1ab58SBarry Smith   PetscInt       nnz          = in_matrix.nnz;
453b2f1ab58SBarry Smith   PetscInt       ncols        = in_matrix.nrows;
454b2f1ab58SBarry Smith   PetscInt       nrows        = in_matrix.ncols;
455b2f1ab58SBarry Smith   PetscInt       i,j,k;
456b2f1ab58SBarry Smith   PetscInt       r_nnz;
457b2f1ab58SBarry Smith   PetscInt       *irow;
4584efc9174SBarry Smith   PetscInt       icol0 = 0;
459b2f1ab58SBarry Smith   PetscScalar    * val;
460b2f1ab58SBarry Smith   PetscErrorCode ierr;
4614e1ff37aSBarry Smith 
462b2f1ab58SBarry Smith   PetscFunctionBegin;
463c328eeadSBarry Smith   /* Copy input values */
464b2f1ab58SBarry Smith   result->nrows        = nrows;
465b2f1ab58SBarry Smith   result->ncols        = ncols;
466b2f1ab58SBarry Smith   result->nnz          = nnz;
467b2f1ab58SBarry Smith   result->col_idx_type = SPBAS_COLUMN_NUMBERS;
468b2f1ab58SBarry Smith   result->block_data   = PETSC_TRUE;
469b2f1ab58SBarry Smith 
470c328eeadSBarry Smith   /* Allocate sparseness pattern */
47171d55d18SBarry Smith   ierr =  spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
472b2f1ab58SBarry Smith 
473c328eeadSBarry Smith   /*  Count the number of nonzeros in each row */
4742205254eSKarl Rupp   for (i = 0; i<nrows; i++) result->row_nnz[i] = 0;
475b2f1ab58SBarry Smith 
4769ccc4a9bSBarry Smith   for (i=0; i<ncols; i++) {
477b2f1ab58SBarry Smith     r_nnz = in_matrix.row_nnz[i];
478b2f1ab58SBarry Smith     irow  = in_matrix.icols[i];
4799ccc4a9bSBarry Smith     if (col_idx_type == SPBAS_COLUMN_NUMBERS)  {
4802205254eSKarl Rupp       for (j=0; j<r_nnz; j++) result->row_nnz[irow[j]]++;
4819ccc4a9bSBarry Smith     } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS)  {
4822205254eSKarl Rupp       for (j=0; j<r_nnz; j++) result->row_nnz[i+irow[j]]++;
4839ccc4a9bSBarry Smith     } else if (col_idx_type == SPBAS_OFFSET_ARRAY) {
484b2f1ab58SBarry Smith       icol0=in_matrix.icol0[i];
4852205254eSKarl Rupp       for (j=0; j<r_nnz; j++) result->row_nnz[icol0+irow[j]]++;
486b2f1ab58SBarry Smith     }
487b2f1ab58SBarry Smith   }
488b2f1ab58SBarry Smith 
489c328eeadSBarry Smith   /* Set the pointers to the data */
490b2f1ab58SBarry Smith   ierr = spbas_allocate_data(result);CHKERRQ(ierr);
491b2f1ab58SBarry Smith 
492c328eeadSBarry Smith   /* Reset the number of nonzeros in each row */
4932205254eSKarl Rupp   for (i = 0; i<nrows; i++) result->row_nnz[i] = 0;
494b2f1ab58SBarry Smith 
495c328eeadSBarry Smith   /* Fill the data arrays */
4969ccc4a9bSBarry Smith   if (in_matrix.values) {
4979ccc4a9bSBarry Smith     for (i=0; i<ncols; i++) {
498b2f1ab58SBarry Smith       r_nnz = in_matrix.row_nnz[i];
499b2f1ab58SBarry Smith       irow  = in_matrix.icols[i];
500b2f1ab58SBarry Smith       val   = in_matrix.values[i];
501b2f1ab58SBarry Smith 
5022205254eSKarl Rupp       if      (col_idx_type == SPBAS_COLUMN_NUMBERS)   icol0 = 0;
5032205254eSKarl Rupp       else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0 = i;
5042205254eSKarl Rupp       else if (col_idx_type == SPBAS_OFFSET_ARRAY)     icol0 = in_matrix.icol0[i];
5059ccc4a9bSBarry Smith       for (j=0; j<r_nnz; j++)  {
506b2f1ab58SBarry Smith         k = icol0 + irow[j];
507b2f1ab58SBarry Smith         result->icols[k][result->row_nnz[k]]  = i;
508b2f1ab58SBarry Smith         result->values[k][result->row_nnz[k]] = val[j];
509b2f1ab58SBarry Smith         result->row_nnz[k]++;
510b2f1ab58SBarry Smith       }
511b2f1ab58SBarry Smith     }
5129ccc4a9bSBarry Smith   } else {
5139ccc4a9bSBarry Smith     for (i=0; i<ncols; i++) {
514b2f1ab58SBarry Smith       r_nnz = in_matrix.row_nnz[i];
515b2f1ab58SBarry Smith       irow  = in_matrix.icols[i];
516b2f1ab58SBarry Smith 
5172205254eSKarl Rupp       if      (col_idx_type == SPBAS_COLUMN_NUMBERS)   icol0=0;
5182205254eSKarl Rupp       else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0=i;
5192205254eSKarl Rupp       else if (col_idx_type == SPBAS_OFFSET_ARRAY)     icol0=in_matrix.icol0[i];
520b2f1ab58SBarry Smith 
5219ccc4a9bSBarry Smith       for (j=0; j<r_nnz; j++) {
522b2f1ab58SBarry Smith         k = icol0 + irow[j];
523b2f1ab58SBarry Smith         result->icols[k][result->row_nnz[k]] = i;
524b2f1ab58SBarry Smith         result->row_nnz[k]++;
525b2f1ab58SBarry Smith       }
526b2f1ab58SBarry Smith     }
527b2f1ab58SBarry Smith   }
528b2f1ab58SBarry Smith   PetscFunctionReturn(0);
529b2f1ab58SBarry Smith }
530b2f1ab58SBarry Smith 
531b2f1ab58SBarry Smith /*
532b2f1ab58SBarry Smith    spbas_mergesort
533b2f1ab58SBarry Smith 
534bb42e86aSJed Brown       mergesort for an array of integers and an array of associated
535b2f1ab58SBarry Smith       reals
536b2f1ab58SBarry Smith 
537b2f1ab58SBarry Smith       on output, icol[0..nnz-1] is increasing;
538b2f1ab58SBarry Smith                   val[0..nnz-1] has undergone the same permutation as icol
539b2f1ab58SBarry Smith 
5400298fd71SBarry Smith       NB: val may be NULL: in that case, only the integers are sorted
541b2f1ab58SBarry Smith 
542b2f1ab58SBarry Smith */
543b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val)
544b2f1ab58SBarry Smith {
545c328eeadSBarry Smith   PetscInt       istep;       /* Chunk-sizes of already sorted parts of arrays */
546c328eeadSBarry Smith   PetscInt       i, i1, i2;   /* Loop counters for (partly) sorted arrays */
547c328eeadSBarry Smith   PetscInt       istart, i1end, i2end; /* start of newly sorted array part, end of both parts */
548c328eeadSBarry Smith   PetscInt       *ialloc;     /* Allocated arrays */
5490298fd71SBarry Smith   PetscScalar    *valloc=NULL;
550c328eeadSBarry Smith   PetscInt       *iswap;      /* auxiliary pointers for swapping */
551b2f1ab58SBarry Smith   PetscScalar    *vswap;
552c328eeadSBarry Smith   PetscInt       *ihlp1;      /* Pointers to new version of arrays, */
5530298fd71SBarry Smith   PetscScalar    *vhlp1=NULL;  /* (arrays under construction) */
554c328eeadSBarry Smith   PetscInt       *ihlp2;      /* Pointers to previous version of arrays, */
5550298fd71SBarry Smith   PetscScalar    *vhlp2=NULL;
556b2f1ab58SBarry Smith   PetscErrorCode ierr;
557b2f1ab58SBarry Smith 
558785e854fSJed Brown   ierr  = PetscMalloc1(nnz,&ialloc);CHKERRQ(ierr);
559b2f1ab58SBarry Smith   ihlp1 = ialloc;
560b2f1ab58SBarry Smith   ihlp2 = icol;
561b2f1ab58SBarry Smith 
5629ccc4a9bSBarry Smith   if (val) {
563785e854fSJed Brown     ierr  = PetscMalloc1(nnz,&valloc);CHKERRQ(ierr);
564b2f1ab58SBarry Smith     vhlp1 = valloc;
565b2f1ab58SBarry Smith     vhlp2 = val;
566b2f1ab58SBarry Smith   }
567b2f1ab58SBarry Smith 
568b2f1ab58SBarry Smith 
569c328eeadSBarry Smith   /* Sorted array chunks are first 1 long, and increase until they are the complete array */
5709ccc4a9bSBarry Smith   for (istep=1; istep<nnz; istep*=2) {
571c328eeadSBarry Smith     /*
572c328eeadSBarry Smith       Combine sorted parts
573c328eeadSBarry Smith           istart:istart+istep-1 and istart+istep-1:istart+2*istep-1
574c328eeadSBarry Smith       of ihlp2 and vhlp2
575c328eeadSBarry Smith 
576c328eeadSBarry Smith       into one sorted part
577c328eeadSBarry Smith           istart:istart+2*istep-1
578c328eeadSBarry Smith       of ihlp1 and vhlp1
579c328eeadSBarry Smith     */
5809ccc4a9bSBarry Smith     for (istart=0; istart<nnz; istart+=2*istep) {
581c328eeadSBarry Smith       /* Set counters and bound array part endings */
5822205254eSKarl Rupp       i1=istart;        i1end = i1+istep;  if (i1end>nnz) i1end=nnz;
5832205254eSKarl Rupp       i2=istart+istep;  i2end = i2+istep;  if (i2end>nnz) i2end=nnz;
584b2f1ab58SBarry Smith 
585c328eeadSBarry Smith       /* Merge the two array parts */
5869ccc4a9bSBarry Smith       if (val) {
5879ccc4a9bSBarry Smith         for (i=istart; i<i2end; i++) {
5889ccc4a9bSBarry Smith           if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) {
589b2f1ab58SBarry Smith             ihlp1[i] = ihlp2[i1];
590b2f1ab58SBarry Smith             vhlp1[i] = vhlp2[i1];
591b2f1ab58SBarry Smith             i1++;
5929ccc4a9bSBarry Smith           } else if (i2<i2end) {
593b2f1ab58SBarry Smith             ihlp1[i] = ihlp2[i2];
594b2f1ab58SBarry Smith             vhlp1[i] = vhlp2[i2];
595b2f1ab58SBarry Smith             i2++;
5969ccc4a9bSBarry Smith           } else {
597b2f1ab58SBarry Smith             ihlp1[i] = ihlp2[i1];
598b2f1ab58SBarry Smith             vhlp1[i] = vhlp2[i1];
599b2f1ab58SBarry Smith             i1++;
600b2f1ab58SBarry Smith           }
601b2f1ab58SBarry Smith         }
6029ccc4a9bSBarry Smith       } else {
6036322f4bdSBarry Smith         for (i=istart; i<i2end; i++) {
6046322f4bdSBarry Smith           if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) {
605b2f1ab58SBarry Smith             ihlp1[i] = ihlp2[i1];
606b2f1ab58SBarry Smith             i1++;
6076322f4bdSBarry Smith           } else if (i2<i2end) {
608b2f1ab58SBarry Smith             ihlp1[i] = ihlp2[i2];
609b2f1ab58SBarry Smith             i2++;
6106322f4bdSBarry Smith           } else {
611b2f1ab58SBarry Smith             ihlp1[i] = ihlp2[i1];
612b2f1ab58SBarry Smith             i1++;
613b2f1ab58SBarry Smith           }
614b2f1ab58SBarry Smith         }
615b2f1ab58SBarry Smith       }
616b2f1ab58SBarry Smith     }
617b2f1ab58SBarry Smith 
618c328eeadSBarry Smith     /* Swap the two array sets */
619b2f1ab58SBarry Smith     iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap;
620b2f1ab58SBarry Smith     vswap = vhlp2; vhlp2 = vhlp1; vhlp1 = vswap;
621b2f1ab58SBarry Smith   }
622b2f1ab58SBarry Smith 
623c328eeadSBarry Smith   /* Copy one more time in case the sorted arrays are the temporary ones */
6246322f4bdSBarry Smith   if (ihlp2 != icol) {
6252205254eSKarl Rupp     for (i=0; i<nnz; i++) icol[i] = ihlp2[i];
6266322f4bdSBarry Smith     if (val) {
6272205254eSKarl Rupp       for (i=0; i<nnz; i++) val[i] = vhlp2[i];
628b2f1ab58SBarry Smith     }
629b2f1ab58SBarry Smith   }
630b2f1ab58SBarry Smith 
631b2f1ab58SBarry Smith   ierr = PetscFree(ialloc);CHKERRQ(ierr);
632b2f1ab58SBarry Smith   if (val) {ierr = PetscFree(valloc);CHKERRQ(ierr);}
633b2f1ab58SBarry Smith   PetscFunctionReturn(0);
634b2f1ab58SBarry Smith }
635b2f1ab58SBarry Smith 
636b2f1ab58SBarry Smith /*
637b2f1ab58SBarry Smith   spbas_apply_reordering_rows:
638b2f1ab58SBarry Smith     apply the given reordering to the rows:  matrix_A = matrix_A(perm,:);
639b2f1ab58SBarry Smith */
640b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_rows(spbas_matrix *matrix_A, const PetscInt *permutation)
641b2f1ab58SBarry Smith {
642b2f1ab58SBarry Smith   PetscInt       i,j,ip;
643b2f1ab58SBarry Smith   PetscInt       nrows=matrix_A->nrows;
644b2f1ab58SBarry Smith   PetscInt       * row_nnz;
645b2f1ab58SBarry Smith   PetscInt       **icols;
646ace3abfcSBarry Smith   PetscBool      do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE;
6470298fd71SBarry Smith   PetscScalar    **vals    = NULL;
648b2f1ab58SBarry Smith   PetscErrorCode ierr;
649b2f1ab58SBarry Smith 
650b2f1ab58SBarry Smith   PetscFunctionBegin;
651e32f2f54SBarry Smith   if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern\n");
652b2f1ab58SBarry Smith 
6536322f4bdSBarry Smith   if (do_values) {
654854ce69bSBarry Smith     ierr = PetscMalloc1(nrows, &vals);CHKERRQ(ierr);
655b2f1ab58SBarry Smith   }
656854ce69bSBarry Smith   ierr = PetscMalloc1(nrows, &row_nnz);CHKERRQ(ierr);
657854ce69bSBarry Smith   ierr = PetscMalloc1(nrows, &icols);CHKERRQ(ierr);
658b2f1ab58SBarry Smith 
6596322f4bdSBarry Smith   for (i=0; i<nrows; i++) {
660b2f1ab58SBarry Smith     ip = permutation[i];
6612205254eSKarl Rupp     if (do_values) vals[i] = matrix_A->values[ip];
662b2f1ab58SBarry Smith     icols[i]   = matrix_A->icols[ip];
663b2f1ab58SBarry Smith     row_nnz[i] = matrix_A->row_nnz[ip];
6642205254eSKarl Rupp     for (j=0; j<row_nnz[i]; j++) icols[i][j] += ip-i;
665b2f1ab58SBarry Smith   }
666b2f1ab58SBarry Smith 
667b2f1ab58SBarry Smith   if (do_values) { ierr = PetscFree(matrix_A->values);CHKERRQ(ierr);}
668b2f1ab58SBarry Smith   ierr = PetscFree(matrix_A->icols);CHKERRQ(ierr);
669b2f1ab58SBarry Smith   ierr = PetscFree(matrix_A->row_nnz);CHKERRQ(ierr);
670b2f1ab58SBarry Smith 
6712205254eSKarl Rupp   if (do_values) matrix_A->values = vals;
672b2f1ab58SBarry Smith   matrix_A->icols   = icols;
673b2f1ab58SBarry Smith   matrix_A->row_nnz = row_nnz;
674b2f1ab58SBarry Smith   PetscFunctionReturn(0);
675b2f1ab58SBarry Smith }
676b2f1ab58SBarry Smith 
677b2f1ab58SBarry Smith 
678b2f1ab58SBarry Smith /*
679b2f1ab58SBarry Smith   spbas_apply_reordering_cols:
680b2f1ab58SBarry Smith     apply the given reordering to the columns:  matrix_A(:,perm) = matrix_A;
681b2f1ab58SBarry Smith */
682b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_cols(spbas_matrix *matrix_A,const PetscInt *permutation)
683b2f1ab58SBarry Smith {
684b2f1ab58SBarry Smith   PetscInt       i,j;
685b2f1ab58SBarry Smith   PetscInt       nrows=matrix_A->nrows;
686b2f1ab58SBarry Smith   PetscInt       row_nnz;
687b2f1ab58SBarry Smith   PetscInt       *icols;
688ace3abfcSBarry Smith   PetscBool      do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE;
6890298fd71SBarry Smith   PetscScalar    *vals     = NULL;
690b2f1ab58SBarry Smith   PetscErrorCode ierr;
691b2f1ab58SBarry Smith 
692b2f1ab58SBarry Smith   PetscFunctionBegin;
693e32f2f54SBarry Smith   if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "must have diagonal offsets in pattern\n");
694b2f1ab58SBarry Smith 
6956322f4bdSBarry Smith   for (i=0; i<nrows; i++) {
696b2f1ab58SBarry Smith     icols   = matrix_A->icols[i];
697b2f1ab58SBarry Smith     row_nnz = matrix_A->row_nnz[i];
6982205254eSKarl Rupp     if (do_values) vals = matrix_A->values[i];
699b2f1ab58SBarry Smith 
7006322f4bdSBarry Smith     for (j=0; j<row_nnz; j++) {
701b2f1ab58SBarry Smith       icols[j] = permutation[i+icols[j]]-i;
702b2f1ab58SBarry Smith     }
703b2f1ab58SBarry Smith     ierr = spbas_mergesort(row_nnz, icols, vals);CHKERRQ(ierr);
704b2f1ab58SBarry Smith   }
705b2f1ab58SBarry Smith   PetscFunctionReturn(0);
706b2f1ab58SBarry Smith }
707b2f1ab58SBarry Smith 
708b2f1ab58SBarry Smith /*
709b2f1ab58SBarry Smith   spbas_apply_reordering:
710b2f1ab58SBarry Smith     apply the given reordering:  matrix_A(perm,perm) = matrix_A;
711b2f1ab58SBarry Smith */
712b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering(spbas_matrix *matrix_A, const PetscInt *permutation, const PetscInt * inv_perm)
713b2f1ab58SBarry Smith {
714b2f1ab58SBarry Smith   PetscErrorCode ierr;
7155fd66863SKarl Rupp 
716b2f1ab58SBarry Smith   PetscFunctionBegin;
717b2f1ab58SBarry Smith   ierr = spbas_apply_reordering_rows(matrix_A, inv_perm);CHKERRQ(ierr);
718b2f1ab58SBarry Smith   ierr = spbas_apply_reordering_cols(matrix_A, permutation);CHKERRQ(ierr);
719b2f1ab58SBarry Smith   PetscFunctionReturn(0);
720b2f1ab58SBarry Smith }
721b2f1ab58SBarry Smith 
722b2f1ab58SBarry Smith PetscErrorCode spbas_pattern_only(PetscInt nrows, PetscInt ncols, PetscInt *ai, PetscInt *aj, spbas_matrix * result)
723b2f1ab58SBarry Smith {
724b2f1ab58SBarry Smith   spbas_matrix   retval;
725b2f1ab58SBarry Smith   PetscInt       i, j, i0, r_nnz;
726b2f1ab58SBarry Smith   PetscErrorCode ierr;
727b2f1ab58SBarry Smith 
728b2f1ab58SBarry Smith   PetscFunctionBegin;
729c328eeadSBarry Smith   /* Copy input values */
730b2f1ab58SBarry Smith   retval.nrows = nrows;
731b2f1ab58SBarry Smith   retval.ncols = ncols;
732b2f1ab58SBarry Smith   retval.nnz   = ai[nrows];
733b2f1ab58SBarry Smith 
734b2f1ab58SBarry Smith   retval.block_data   = PETSC_TRUE;
735b2f1ab58SBarry Smith   retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS;
736b2f1ab58SBarry Smith 
737c328eeadSBarry Smith   /* Allocate output matrix */
738b2f1ab58SBarry Smith   ierr = spbas_allocate_pattern(&retval, PETSC_FALSE);CHKERRQ(ierr);
7392205254eSKarl Rupp   for (i=0; i<nrows; i++) retval.row_nnz[i] = ai[i+1]-ai[i];
740b2f1ab58SBarry Smith   ierr = spbas_allocate_data(&retval);CHKERRQ(ierr);
741c328eeadSBarry Smith   /* Copy the structure */
7426322f4bdSBarry Smith   for (i = 0; i<retval.nrows; i++)  {
743b2f1ab58SBarry Smith     i0    = ai[i];
744b2f1ab58SBarry Smith     r_nnz = ai[i+1]-i0;
745b2f1ab58SBarry Smith 
7466322f4bdSBarry Smith     for (j=0; j<r_nnz; j++) {
747b2f1ab58SBarry Smith       retval.icols[i][j] = aj[i0+j]-i;
748b2f1ab58SBarry Smith     }
749b2f1ab58SBarry Smith   }
750b2f1ab58SBarry Smith   *result = retval;
751b2f1ab58SBarry Smith   PetscFunctionReturn(0);
752b2f1ab58SBarry Smith }
753b2f1ab58SBarry Smith 
754b2f1ab58SBarry Smith 
755b2f1ab58SBarry Smith /*
756b2f1ab58SBarry Smith    spbas_mark_row_power:
757b2f1ab58SBarry Smith       Mark the columns in row 'row' which are nonzero in
758b2f1ab58SBarry Smith           matrix^2log(marker).
759b2f1ab58SBarry Smith */
760be332245SKarl Rupp PetscErrorCode spbas_mark_row_power(PetscInt *iwork,             /* marker-vector */
761c328eeadSBarry Smith                                     PetscInt row,                /* row for which the columns are marked */
762c328eeadSBarry Smith                                     spbas_matrix * in_matrix,    /* matrix for which the power is being  calculated */
763c328eeadSBarry Smith                                     PetscInt marker,             /* marker-value: 2^power */
764c328eeadSBarry Smith                                     PetscInt minmrk,             /* lower bound for marked points */
765c328eeadSBarry Smith                                     PetscInt maxmrk)             /* upper bound for marked points */
766b2f1ab58SBarry Smith {
767b2f1ab58SBarry Smith   PetscErrorCode ierr;
768b2f1ab58SBarry Smith   PetscInt       i,j, nnz;
769b2f1ab58SBarry Smith 
770b2f1ab58SBarry Smith   PetscFunctionBegin;
771b2f1ab58SBarry Smith   nnz = in_matrix->row_nnz[row];
772b2f1ab58SBarry Smith 
773c328eeadSBarry Smith   /* For higher powers, call this function recursively */
7746322f4bdSBarry Smith   if (marker>1) {
7756322f4bdSBarry Smith     for (i=0; i<nnz; i++) {
776b2f1ab58SBarry Smith       j = row + in_matrix->icols[row][i];
7776322f4bdSBarry Smith       if (minmrk<=j && j<maxmrk && iwork[j] < marker) {
77871d55d18SBarry Smith         ierr      = spbas_mark_row_power(iwork, row + in_matrix->icols[row][i],in_matrix, marker/2,minmrk,maxmrk);CHKERRQ(ierr);
779b2f1ab58SBarry Smith         iwork[j] |= marker;
780b2f1ab58SBarry Smith       }
781b2f1ab58SBarry Smith     }
7826322f4bdSBarry Smith   } else {
783c328eeadSBarry Smith     /*  Mark the columns reached */
7846322f4bdSBarry Smith     for (i=0; i<nnz; i++)  {
785b2f1ab58SBarry Smith       j = row + in_matrix->icols[row][i];
7862205254eSKarl Rupp       if (minmrk<=j && j<maxmrk) iwork[j] |= 1;
787b2f1ab58SBarry Smith     }
788b2f1ab58SBarry Smith   }
789b2f1ab58SBarry Smith   PetscFunctionReturn(0);
790b2f1ab58SBarry Smith }
791b2f1ab58SBarry Smith 
792b2f1ab58SBarry Smith 
793b2f1ab58SBarry Smith /*
794b2f1ab58SBarry Smith    spbas_power
795b2f1ab58SBarry Smith       Calculate sparseness patterns for incomplete Cholesky decompositions
796b2f1ab58SBarry Smith       of a given order: (almost) all nonzeros of the matrix^(order+1) which
797b2f1ab58SBarry Smith       are inside the band width are found and stored in the output sparseness
798b2f1ab58SBarry Smith       pattern.
799b2f1ab58SBarry Smith */
800b2f1ab58SBarry Smith PetscErrorCode spbas_power(spbas_matrix in_matrix,PetscInt power, spbas_matrix * result)
801b2f1ab58SBarry Smith {
802b2f1ab58SBarry Smith   spbas_matrix   retval;
803b2f1ab58SBarry Smith   PetscInt       nrows = in_matrix.nrows;
804b2f1ab58SBarry Smith   PetscInt       ncols = in_matrix.ncols;
805b2f1ab58SBarry Smith   PetscInt       i, j, kend;
806b2f1ab58SBarry Smith   PetscInt       nnz, inz;
807b2f1ab58SBarry Smith   PetscInt       *iwork;
808b2f1ab58SBarry Smith   PetscInt       marker;
809b2f1ab58SBarry Smith   PetscInt       maxmrk=0;
810b2f1ab58SBarry Smith   PetscErrorCode ierr;
811b2f1ab58SBarry Smith 
812b2f1ab58SBarry Smith   PetscFunctionBegin;
813e32f2f54SBarry Smith   if (in_matrix.col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern\n");
814e32f2f54SBarry Smith   if (ncols != nrows) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Dimension error\n");
815e32f2f54SBarry Smith   if (in_matrix.values) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input array must be sparseness pattern (no values)");
816e32f2f54SBarry Smith   if (power<=0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Power must be 1 or up");
817b2f1ab58SBarry Smith 
818c328eeadSBarry Smith   /* Copy input values*/
819b2f1ab58SBarry Smith   retval.nrows        = ncols;
820b2f1ab58SBarry Smith   retval.ncols        = nrows;
821b2f1ab58SBarry Smith   retval.nnz          = 0;
822b2f1ab58SBarry Smith   retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS;
823b2f1ab58SBarry Smith   retval.block_data   = PETSC_FALSE;
824b2f1ab58SBarry Smith 
825c328eeadSBarry Smith   /* Allocate sparseness pattern */
82671d55d18SBarry Smith   ierr =  spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
827b2f1ab58SBarry Smith 
828*580bdb30SBarry Smith   /* Allocate marker array: note sure the max needed so use the max of the two */
829*580bdb30SBarry Smith   ierr = PetscCalloc1(PetscMax(ncols,nrows), &iwork);CHKERRQ(ierr);
830b2f1ab58SBarry Smith 
831c328eeadSBarry Smith   /* Calculate marker values */
8322205254eSKarl Rupp   marker = 1; for (i=1; i<power; i++) marker*=2;
833b2f1ab58SBarry Smith 
8346322f4bdSBarry Smith   for (i=0; i<nrows; i++)  {
835c328eeadSBarry Smith     /* Calculate the pattern for each row */
836b2f1ab58SBarry Smith 
837b2f1ab58SBarry Smith     nnz  = in_matrix.row_nnz[i];
838b2f1ab58SBarry Smith     kend = i+in_matrix.icols[i][nnz-1];
8392205254eSKarl Rupp     if (maxmrk<=kend) maxmrk=kend+1;
84071d55d18SBarry Smith     ierr = spbas_mark_row_power(iwork, i, &in_matrix, marker, i, maxmrk);CHKERRQ(ierr);
841b2f1ab58SBarry Smith 
842c328eeadSBarry Smith     /* Count the columns*/
843b2f1ab58SBarry Smith     nnz = 0;
8442205254eSKarl Rupp     for (j=i; j<maxmrk; j++) nnz+= (iwork[j]!=0);
845b2f1ab58SBarry Smith 
846c328eeadSBarry Smith     /* Allocate the column indices */
847b2f1ab58SBarry Smith     retval.row_nnz[i] = nnz;
848785e854fSJed Brown     ierr = PetscMalloc1(nnz,&retval.icols[i]);CHKERRQ(ierr);
849b2f1ab58SBarry Smith 
850c328eeadSBarry Smith     /* Administrate the column indices */
851b2f1ab58SBarry Smith     inz = 0;
8526322f4bdSBarry Smith     for (j=i; j<maxmrk; j++) {
8536322f4bdSBarry Smith       if (iwork[j]) {
854b2f1ab58SBarry Smith         retval.icols[i][inz] = j-i;
855b2f1ab58SBarry Smith         inz++;
856b2f1ab58SBarry Smith         iwork[j] = 0;
857b2f1ab58SBarry Smith       }
858b2f1ab58SBarry Smith     }
859b2f1ab58SBarry Smith     retval.nnz += nnz;
860b2f1ab58SBarry Smith   };
861b2f1ab58SBarry Smith   ierr    = PetscFree(iwork);CHKERRQ(ierr);
862b2f1ab58SBarry Smith   *result = retval;
863b2f1ab58SBarry Smith   PetscFunctionReturn(0);
864b2f1ab58SBarry Smith }
865b2f1ab58SBarry Smith 
866b2f1ab58SBarry Smith 
867b2f1ab58SBarry Smith 
868b2f1ab58SBarry Smith /*
869b2f1ab58SBarry Smith    spbas_keep_upper:
870b2f1ab58SBarry Smith       remove the lower part of the matrix: keep the upper part
871b2f1ab58SBarry Smith */
872b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper(spbas_matrix * inout_matrix)
873b2f1ab58SBarry Smith {
874b2f1ab58SBarry Smith   PetscInt i, j;
875b2f1ab58SBarry Smith   PetscInt jstart;
876b2f1ab58SBarry Smith 
877b2f1ab58SBarry Smith   PetscFunctionBegin;
878e32f2f54SBarry Smith   if (inout_matrix->block_data) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Not yet for block data matrices\n");
8796322f4bdSBarry Smith   for (i=0; i<inout_matrix->nrows; i++)  {
8806322f4bdSBarry Smith     for (jstart=0; (jstart<inout_matrix->row_nnz[i]) && (inout_matrix->icols[i][jstart]<0); jstart++) {}
8816322f4bdSBarry Smith     if (jstart>0) {
8826322f4bdSBarry Smith       for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) {
8836322f4bdSBarry Smith         inout_matrix->icols[i][j] = inout_matrix->icols[i][j+jstart];
884b2f1ab58SBarry Smith       }
885b2f1ab58SBarry Smith 
8866322f4bdSBarry Smith       if (inout_matrix->values) {
8876322f4bdSBarry Smith         for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) {
8886322f4bdSBarry Smith           inout_matrix->values[i][j] = inout_matrix->values[i][j+jstart];
889b2f1ab58SBarry Smith         }
890b2f1ab58SBarry Smith       }
891b2f1ab58SBarry Smith 
892b2f1ab58SBarry Smith       inout_matrix->row_nnz[i] -= jstart;
893b2f1ab58SBarry Smith 
8946322f4bdSBarry Smith       inout_matrix->icols[i] = (PetscInt*) realloc((void*) inout_matrix->icols[i], inout_matrix->row_nnz[i]*sizeof(PetscInt));
895b2f1ab58SBarry Smith 
8966322f4bdSBarry Smith       if (inout_matrix->values) {
8976322f4bdSBarry Smith         inout_matrix->values[i] = (PetscScalar*) realloc((void*) inout_matrix->values[i], inout_matrix->row_nnz[i]*sizeof(PetscScalar));
898b2f1ab58SBarry Smith       }
899b2f1ab58SBarry Smith       inout_matrix->nnz -= jstart;
900b2f1ab58SBarry Smith     }
901b2f1ab58SBarry Smith   }
902b2f1ab58SBarry Smith   PetscFunctionReturn(0);
903b2f1ab58SBarry Smith }
904b2f1ab58SBarry Smith 
905