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 spbas_compress_pattern: 248b2f1ab58SBarry Smith calculate a compressed sparseness pattern for a sparseness pattern 249b2f1ab58SBarry Smith given in compressed row storage. The compressed sparseness pattern may 250b2f1ab58SBarry Smith require (much) less memory. 251b2f1ab58SBarry Smith */ 2524e1ff37aSBarry 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) 253b2f1ab58SBarry Smith { 254b2f1ab58SBarry Smith PetscInt nnz = irow_in[nrows]; 2553dfa2556SBarry Smith size_t mem_orig = (nrows + nnz) * sizeof(PetscInt); 2563dfa2556SBarry Smith size_t mem_compressed; 257b2f1ab58SBarry Smith PetscErrorCode ierr; 258b2f1ab58SBarry Smith PetscInt *isort; 259b2f1ab58SBarry Smith PetscInt *icols; 260b2f1ab58SBarry Smith PetscInt row_nnz; 261b2f1ab58SBarry Smith PetscInt *ipoint; 262ace3abfcSBarry Smith PetscBool *used; 263b2f1ab58SBarry Smith PetscInt ptr; 264b2f1ab58SBarry Smith PetscInt i,j; 265ace3abfcSBarry Smith const PetscBool no_values = PETSC_FALSE; 266b2f1ab58SBarry Smith 267b2f1ab58SBarry Smith PetscFunctionBegin; 268c328eeadSBarry Smith /* Allocate the structure of the new matrix */ 269b2f1ab58SBarry Smith B->nrows = nrows; 270b2f1ab58SBarry Smith B->ncols = ncols; 271b2f1ab58SBarry Smith B->nnz = nnz; 272b2f1ab58SBarry Smith B->col_idx_type = col_idx_type; 273b2f1ab58SBarry Smith B->block_data = PETSC_TRUE; 2742205254eSKarl Rupp 275b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(B, no_values);CHKERRQ(ierr); 276b2f1ab58SBarry Smith 277c328eeadSBarry Smith /* When using an offset array, set it */ 2784e1ff37aSBarry Smith if (col_idx_type==SPBAS_OFFSET_ARRAY) { 2792205254eSKarl Rupp for (i=0; i<nrows; i++) B->icol0[i] = icol_in[irow_in[i]]; 280b2f1ab58SBarry Smith } 281b2f1ab58SBarry Smith 282c328eeadSBarry Smith /* Allocate the ordering for the rows */ 283785e854fSJed Brown ierr = PetscMalloc1(nrows,&isort);CHKERRQ(ierr); 284785e854fSJed Brown ierr = PetscMalloc1(nrows,&ipoint);CHKERRQ(ierr); 285580bdb30SBarry Smith ierr = PetscCalloc1(nrows,&used);CHKERRQ(ierr); 286b2f1ab58SBarry Smith 2874e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 288b2f1ab58SBarry Smith B->row_nnz[i] = irow_in[i+1]-irow_in[i]; 289b2f1ab58SBarry Smith isort[i] = i; 290b2f1ab58SBarry Smith ipoint[i] = i; 291b2f1ab58SBarry Smith } 292b2f1ab58SBarry Smith 293c328eeadSBarry Smith /* Sort the rows so that identical columns will be next to each other */ 294b2f1ab58SBarry Smith ierr = spbas_mergesort_icols(nrows, irow_in, icol_in, col_idx_type, isort);CHKERRQ(ierr); 2950298fd71SBarry Smith ierr = PetscInfo(NULL,"Rows have been sorted for patterns\n");CHKERRQ(ierr); 296b2f1ab58SBarry Smith 297c328eeadSBarry Smith /* Replace identical rows with the first one in the list */ 2984e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 2994e1ff37aSBarry Smith if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) { 300b2f1ab58SBarry Smith ipoint[isort[i]] = ipoint[isort[i-1]]; 301b2f1ab58SBarry Smith } 302b2f1ab58SBarry Smith } 303b2f1ab58SBarry Smith 304c328eeadSBarry Smith /* Collect the rows which are used*/ 3052205254eSKarl Rupp for (i=0; i<nrows; i++) used[ipoint[i]] = PETSC_TRUE; 306b2f1ab58SBarry Smith 307c328eeadSBarry Smith /* Calculate needed memory */ 308b2f1ab58SBarry Smith B->n_alloc_icol = 0; 3094e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 3102205254eSKarl Rupp if (used[i]) B->n_alloc_icol += B->row_nnz[i]; 311b2f1ab58SBarry Smith } 312785e854fSJed Brown ierr = PetscMalloc1(B->n_alloc_icol,&B->alloc_icol);CHKERRQ(ierr); 313b2f1ab58SBarry Smith 314c328eeadSBarry Smith /* Fill in the diagonal offsets for the rows which store their own data */ 315b2f1ab58SBarry Smith ptr = 0; 3164e1ff37aSBarry Smith for (i=0; i<B->nrows; i++) { 3174e1ff37aSBarry Smith if (used[i]) { 318b2f1ab58SBarry Smith B->icols[i] = &B->alloc_icol[ptr]; 319b2f1ab58SBarry Smith icols = &icol_in[irow_in[i]]; 320b2f1ab58SBarry Smith row_nnz = B->row_nnz[i]; 3214e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 3224e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 323b2f1ab58SBarry Smith B->icols[i][j] = icols[j]; 324b2f1ab58SBarry Smith } 3254e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 3264e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 327b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-i; 328b2f1ab58SBarry Smith } 3294e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 3304e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 331b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-icols[0]; 332b2f1ab58SBarry Smith } 333b2f1ab58SBarry Smith } 334b2f1ab58SBarry Smith ptr += B->row_nnz[i]; 335b2f1ab58SBarry Smith } 336b2f1ab58SBarry Smith } 337b2f1ab58SBarry Smith 338c328eeadSBarry Smith /* Point to the right places for all data */ 3394e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 340b2f1ab58SBarry Smith B->icols[i] = B->icols[ipoint[i]]; 341b2f1ab58SBarry Smith } 3420298fd71SBarry Smith ierr = PetscInfo(NULL,"Row patterns have been compressed\n");CHKERRQ(ierr); 34357622a8eSBarry Smith ierr = PetscInfo1(NULL," (%g nonzeros per row)\n", (double) ((PetscReal) nnz / (PetscReal) nrows));CHKERRQ(ierr); 344b2f1ab58SBarry Smith 345b2f1ab58SBarry Smith ierr=PetscFree(isort);CHKERRQ(ierr); 346b2f1ab58SBarry Smith ierr=PetscFree(used);CHKERRQ(ierr); 347b2f1ab58SBarry Smith ierr=PetscFree(ipoint);CHKERRQ(ierr); 348b2f1ab58SBarry Smith 349b2f1ab58SBarry Smith mem_compressed = spbas_memory_requirement(*B); 3504e1ff37aSBarry Smith *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig; 351b2f1ab58SBarry Smith PetscFunctionReturn(0); 352b2f1ab58SBarry Smith } 353b2f1ab58SBarry Smith 354b2f1ab58SBarry Smith /* 355b2f1ab58SBarry Smith spbas_incomplete_cholesky 356b2f1ab58SBarry Smith Incomplete Cholesky decomposition 357b2f1ab58SBarry Smith */ 358c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas_cholesky.h> 359b2f1ab58SBarry Smith 360b2f1ab58SBarry Smith /* 361b2f1ab58SBarry Smith spbas_delete : de-allocate the arrays owned by this matrix 362b2f1ab58SBarry Smith */ 363b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix matrix) 364b2f1ab58SBarry Smith { 365b2f1ab58SBarry Smith PetscInt i; 366b2f1ab58SBarry Smith PetscErrorCode ierr; 3679ccc4a9bSBarry Smith 368b2f1ab58SBarry Smith PetscFunctionBegin; 3699ccc4a9bSBarry Smith if (matrix.block_data) { 370cb9801acSJed Brown ierr=PetscFree(matrix.alloc_icol);CHKERRQ(ierr); 371b2f1ab58SBarry Smith if (matrix.values) {ierr=PetscFree(matrix.alloc_val);CHKERRQ(ierr);} 3729ccc4a9bSBarry Smith } else { 3739ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.icols[i]);CHKERRQ(ierr);} 374b2f1ab58SBarry Smith ierr = PetscFree(matrix.icols);CHKERRQ(ierr); 3759ccc4a9bSBarry Smith if (matrix.values) { 3769ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.values[i]);CHKERRQ(ierr);} 377b2f1ab58SBarry Smith } 378b2f1ab58SBarry Smith } 379b2f1ab58SBarry Smith 380b2f1ab58SBarry Smith ierr=PetscFree(matrix.row_nnz);CHKERRQ(ierr); 381b2f1ab58SBarry Smith ierr=PetscFree(matrix.icols);CHKERRQ(ierr); 3829ccc4a9bSBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) {ierr=PetscFree(matrix.icol0);CHKERRQ(ierr);} 383c31cb41cSBarry Smith ierr=PetscFree(matrix.values);CHKERRQ(ierr); 384b2f1ab58SBarry Smith PetscFunctionReturn(0); 385b2f1ab58SBarry Smith } 386b2f1ab58SBarry Smith 387b2f1ab58SBarry Smith /* 388b2f1ab58SBarry Smith spbas_matrix_to_crs: 389b2f1ab58SBarry Smith Convert an spbas_matrix to compessed row storage 390b2f1ab58SBarry Smith */ 391b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out) 392b2f1ab58SBarry Smith { 393b2f1ab58SBarry Smith PetscInt nrows = matrix_A.nrows; 394b2f1ab58SBarry Smith PetscInt nnz = matrix_A.nnz; 395b2f1ab58SBarry Smith PetscInt i,j,r_nnz,i0; 396b2f1ab58SBarry Smith PetscInt *irow; 397b2f1ab58SBarry Smith PetscInt *icol; 398b2f1ab58SBarry Smith PetscInt *icol_A; 399b2f1ab58SBarry Smith MatScalar *val; 400b2f1ab58SBarry Smith PetscScalar *val_A; 401b2f1ab58SBarry Smith PetscInt col_idx_type = matrix_A.col_idx_type; 402ace3abfcSBarry Smith PetscBool do_values = matrix_A.values ? PETSC_TRUE : PETSC_FALSE; 403b2f1ab58SBarry Smith PetscErrorCode ierr; 404b2f1ab58SBarry Smith 405b2f1ab58SBarry Smith PetscFunctionBegin; 406854ce69bSBarry Smith ierr = PetscMalloc1(nrows+1, &irow);CHKERRQ(ierr); 407854ce69bSBarry Smith ierr = PetscMalloc1(nnz, &icol);CHKERRQ(ierr); 4089ccc4a9bSBarry Smith *icol_out = icol; 4099ccc4a9bSBarry Smith *irow_out = irow; 4109ccc4a9bSBarry Smith if (do_values) { 411854ce69bSBarry Smith ierr = PetscMalloc1(nnz, &val);CHKERRQ(ierr); 412b2f1ab58SBarry Smith *val_out = val; *icol_out = icol; *irow_out=irow; 413b2f1ab58SBarry Smith } 414b2f1ab58SBarry Smith 415b2f1ab58SBarry Smith irow[0]=0; 4169ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 417b2f1ab58SBarry Smith r_nnz = matrix_A.row_nnz[i]; 418b2f1ab58SBarry Smith i0 = irow[i]; 419b2f1ab58SBarry Smith irow[i+1] = i0 + r_nnz; 420b2f1ab58SBarry Smith icol_A = matrix_A.icols[i]; 421b2f1ab58SBarry Smith 4229ccc4a9bSBarry Smith if (do_values) { 423b2f1ab58SBarry Smith val_A = matrix_A.values[i]; 4249ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 425b2f1ab58SBarry Smith icol[i0+j] = icol_A[j]; 426b2f1ab58SBarry Smith val[i0+j] = val_A[j]; 427b2f1ab58SBarry Smith } 4289ccc4a9bSBarry Smith } else { 4292205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] = icol_A[j]; 430b2f1ab58SBarry Smith } 431b2f1ab58SBarry Smith 4329ccc4a9bSBarry Smith if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 4332205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] += i; 4342205254eSKarl Rupp } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 435b2f1ab58SBarry Smith i0 = matrix_A.icol0[i]; 4362205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] += i0; 437b2f1ab58SBarry Smith } 438b2f1ab58SBarry Smith } 439b2f1ab58SBarry Smith PetscFunctionReturn(0); 440b2f1ab58SBarry Smith } 441b2f1ab58SBarry Smith 442b2f1ab58SBarry Smith /* 443b2f1ab58SBarry Smith spbas_transpose 444b2f1ab58SBarry Smith return the transpose of a matrix 445b2f1ab58SBarry Smith */ 446b2f1ab58SBarry Smith PetscErrorCode spbas_transpose(spbas_matrix in_matrix, spbas_matrix * result) 447b2f1ab58SBarry Smith { 448b2f1ab58SBarry Smith PetscInt col_idx_type = in_matrix.col_idx_type; 449b2f1ab58SBarry Smith PetscInt nnz = in_matrix.nnz; 450b2f1ab58SBarry Smith PetscInt ncols = in_matrix.nrows; 451b2f1ab58SBarry Smith PetscInt nrows = in_matrix.ncols; 452b2f1ab58SBarry Smith PetscInt i,j,k; 453b2f1ab58SBarry Smith PetscInt r_nnz; 454b2f1ab58SBarry Smith PetscInt *irow; 4554efc9174SBarry Smith PetscInt icol0 = 0; 456b2f1ab58SBarry Smith PetscScalar * val; 457b2f1ab58SBarry Smith PetscErrorCode ierr; 4584e1ff37aSBarry Smith 459b2f1ab58SBarry Smith PetscFunctionBegin; 460c328eeadSBarry Smith /* Copy input values */ 461b2f1ab58SBarry Smith result->nrows = nrows; 462b2f1ab58SBarry Smith result->ncols = ncols; 463b2f1ab58SBarry Smith result->nnz = nnz; 464b2f1ab58SBarry Smith result->col_idx_type = SPBAS_COLUMN_NUMBERS; 465b2f1ab58SBarry Smith result->block_data = PETSC_TRUE; 466b2f1ab58SBarry Smith 467c328eeadSBarry Smith /* Allocate sparseness pattern */ 46871d55d18SBarry Smith ierr = spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 469b2f1ab58SBarry Smith 470c328eeadSBarry Smith /* Count the number of nonzeros in each row */ 4712205254eSKarl Rupp for (i = 0; i<nrows; i++) result->row_nnz[i] = 0; 472b2f1ab58SBarry Smith 4739ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 474b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 475b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 4769ccc4a9bSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 4772205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[irow[j]]++; 4789ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 4792205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[i+irow[j]]++; 4809ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 481b2f1ab58SBarry Smith icol0=in_matrix.icol0[i]; 4822205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[icol0+irow[j]]++; 483b2f1ab58SBarry Smith } 484b2f1ab58SBarry Smith } 485b2f1ab58SBarry Smith 486c328eeadSBarry Smith /* Set the pointers to the data */ 487b2f1ab58SBarry Smith ierr = spbas_allocate_data(result);CHKERRQ(ierr); 488b2f1ab58SBarry Smith 489c328eeadSBarry Smith /* Reset the number of nonzeros in each row */ 4902205254eSKarl Rupp for (i = 0; i<nrows; i++) result->row_nnz[i] = 0; 491b2f1ab58SBarry Smith 492c328eeadSBarry Smith /* Fill the data arrays */ 4939ccc4a9bSBarry Smith if (in_matrix.values) { 4949ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 495b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 496b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 497b2f1ab58SBarry Smith val = in_matrix.values[i]; 498b2f1ab58SBarry Smith 4992205254eSKarl Rupp if (col_idx_type == SPBAS_COLUMN_NUMBERS) icol0 = 0; 5002205254eSKarl Rupp else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0 = i; 5012205254eSKarl Rupp else if (col_idx_type == SPBAS_OFFSET_ARRAY) icol0 = in_matrix.icol0[i]; 5029ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 503b2f1ab58SBarry Smith k = icol0 + irow[j]; 504b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 505b2f1ab58SBarry Smith result->values[k][result->row_nnz[k]] = val[j]; 506b2f1ab58SBarry Smith result->row_nnz[k]++; 507b2f1ab58SBarry Smith } 508b2f1ab58SBarry Smith } 5099ccc4a9bSBarry Smith } else { 5109ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 511b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 512b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 513b2f1ab58SBarry Smith 5142205254eSKarl Rupp if (col_idx_type == SPBAS_COLUMN_NUMBERS) icol0=0; 5152205254eSKarl Rupp else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0=i; 5162205254eSKarl Rupp else if (col_idx_type == SPBAS_OFFSET_ARRAY) icol0=in_matrix.icol0[i]; 517b2f1ab58SBarry Smith 5189ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 519b2f1ab58SBarry Smith k = icol0 + irow[j]; 520b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 521b2f1ab58SBarry Smith result->row_nnz[k]++; 522b2f1ab58SBarry Smith } 523b2f1ab58SBarry Smith } 524b2f1ab58SBarry Smith } 525b2f1ab58SBarry Smith PetscFunctionReturn(0); 526b2f1ab58SBarry Smith } 527b2f1ab58SBarry Smith 528b2f1ab58SBarry Smith /* 529b2f1ab58SBarry Smith spbas_mergesort 530b2f1ab58SBarry Smith 531bb42e86aSJed Brown mergesort for an array of integers and an array of associated 532b2f1ab58SBarry Smith reals 533b2f1ab58SBarry Smith 534b2f1ab58SBarry Smith on output, icol[0..nnz-1] is increasing; 535b2f1ab58SBarry Smith val[0..nnz-1] has undergone the same permutation as icol 536b2f1ab58SBarry Smith 5370298fd71SBarry Smith NB: val may be NULL: in that case, only the integers are sorted 538b2f1ab58SBarry Smith 539b2f1ab58SBarry Smith */ 540b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val) 541b2f1ab58SBarry Smith { 542c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 543c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 544c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 545c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 5460298fd71SBarry Smith PetscScalar *valloc=NULL; 547c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 548b2f1ab58SBarry Smith PetscScalar *vswap; 549c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 5500298fd71SBarry Smith PetscScalar *vhlp1=NULL; /* (arrays under construction) */ 551c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 5520298fd71SBarry Smith PetscScalar *vhlp2=NULL; 553b2f1ab58SBarry Smith PetscErrorCode ierr; 554b2f1ab58SBarry Smith 555*362febeeSStefano Zampini PetscFunctionBegin; 556785e854fSJed Brown ierr = PetscMalloc1(nnz,&ialloc);CHKERRQ(ierr); 557b2f1ab58SBarry Smith ihlp1 = ialloc; 558b2f1ab58SBarry Smith ihlp2 = icol; 559b2f1ab58SBarry Smith 5609ccc4a9bSBarry Smith if (val) { 561785e854fSJed Brown ierr = PetscMalloc1(nnz,&valloc);CHKERRQ(ierr); 562b2f1ab58SBarry Smith vhlp1 = valloc; 563b2f1ab58SBarry Smith vhlp2 = val; 564b2f1ab58SBarry Smith } 565b2f1ab58SBarry Smith 566c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 5679ccc4a9bSBarry Smith for (istep=1; istep<nnz; istep*=2) { 568c328eeadSBarry Smith /* 569c328eeadSBarry Smith Combine sorted parts 570c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 571c328eeadSBarry Smith of ihlp2 and vhlp2 572c328eeadSBarry Smith 573c328eeadSBarry Smith into one sorted part 574c328eeadSBarry Smith istart:istart+2*istep-1 575c328eeadSBarry Smith of ihlp1 and vhlp1 576c328eeadSBarry Smith */ 5779ccc4a9bSBarry Smith for (istart=0; istart<nnz; istart+=2*istep) { 578c328eeadSBarry Smith /* Set counters and bound array part endings */ 5792205254eSKarl Rupp i1=istart; i1end = i1+istep; if (i1end>nnz) i1end=nnz; 5802205254eSKarl Rupp i2=istart+istep; i2end = i2+istep; if (i2end>nnz) i2end=nnz; 581b2f1ab58SBarry Smith 582c328eeadSBarry Smith /* Merge the two array parts */ 5839ccc4a9bSBarry Smith if (val) { 5849ccc4a9bSBarry Smith for (i=istart; i<i2end; i++) { 5859ccc4a9bSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 586b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 587b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 588b2f1ab58SBarry Smith i1++; 5899ccc4a9bSBarry Smith } else if (i2<i2end) { 590b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 591b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i2]; 592b2f1ab58SBarry Smith i2++; 5939ccc4a9bSBarry Smith } else { 594b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 595b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 596b2f1ab58SBarry Smith i1++; 597b2f1ab58SBarry Smith } 598b2f1ab58SBarry Smith } 5999ccc4a9bSBarry Smith } else { 6006322f4bdSBarry Smith for (i=istart; i<i2end; i++) { 6016322f4bdSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 602b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 603b2f1ab58SBarry Smith i1++; 6046322f4bdSBarry Smith } else if (i2<i2end) { 605b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 606b2f1ab58SBarry Smith i2++; 6076322f4bdSBarry Smith } else { 608b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 609b2f1ab58SBarry Smith i1++; 610b2f1ab58SBarry Smith } 611b2f1ab58SBarry Smith } 612b2f1ab58SBarry Smith } 613b2f1ab58SBarry Smith } 614b2f1ab58SBarry Smith 615c328eeadSBarry Smith /* Swap the two array sets */ 616b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 617b2f1ab58SBarry Smith vswap = vhlp2; vhlp2 = vhlp1; vhlp1 = vswap; 618b2f1ab58SBarry Smith } 619b2f1ab58SBarry Smith 620c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 6216322f4bdSBarry Smith if (ihlp2 != icol) { 6222205254eSKarl Rupp for (i=0; i<nnz; i++) icol[i] = ihlp2[i]; 6236322f4bdSBarry Smith if (val) { 6242205254eSKarl Rupp for (i=0; i<nnz; i++) val[i] = vhlp2[i]; 625b2f1ab58SBarry Smith } 626b2f1ab58SBarry Smith } 627b2f1ab58SBarry Smith 628b2f1ab58SBarry Smith ierr = PetscFree(ialloc);CHKERRQ(ierr); 629b2f1ab58SBarry Smith if (val) {ierr = PetscFree(valloc);CHKERRQ(ierr);} 630b2f1ab58SBarry Smith PetscFunctionReturn(0); 631b2f1ab58SBarry Smith } 632b2f1ab58SBarry Smith 633b2f1ab58SBarry Smith /* 634b2f1ab58SBarry Smith spbas_apply_reordering_rows: 635b2f1ab58SBarry Smith apply the given reordering to the rows: matrix_A = matrix_A(perm,:); 636b2f1ab58SBarry Smith */ 637b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_rows(spbas_matrix *matrix_A, const PetscInt *permutation) 638b2f1ab58SBarry Smith { 639b2f1ab58SBarry Smith PetscInt i,j,ip; 640b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 641b2f1ab58SBarry Smith PetscInt * row_nnz; 642b2f1ab58SBarry Smith PetscInt **icols; 643ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 6440298fd71SBarry Smith PetscScalar **vals = NULL; 645b2f1ab58SBarry Smith PetscErrorCode ierr; 646b2f1ab58SBarry Smith 647b2f1ab58SBarry Smith PetscFunctionBegin; 648e32f2f54SBarry 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"); 649b2f1ab58SBarry Smith 6506322f4bdSBarry Smith if (do_values) { 651854ce69bSBarry Smith ierr = PetscMalloc1(nrows, &vals);CHKERRQ(ierr); 652b2f1ab58SBarry Smith } 653854ce69bSBarry Smith ierr = PetscMalloc1(nrows, &row_nnz);CHKERRQ(ierr); 654854ce69bSBarry Smith ierr = PetscMalloc1(nrows, &icols);CHKERRQ(ierr); 655b2f1ab58SBarry Smith 6566322f4bdSBarry Smith for (i=0; i<nrows; i++) { 657b2f1ab58SBarry Smith ip = permutation[i]; 6582205254eSKarl Rupp if (do_values) vals[i] = matrix_A->values[ip]; 659b2f1ab58SBarry Smith icols[i] = matrix_A->icols[ip]; 660b2f1ab58SBarry Smith row_nnz[i] = matrix_A->row_nnz[ip]; 6612205254eSKarl Rupp for (j=0; j<row_nnz[i]; j++) icols[i][j] += ip-i; 662b2f1ab58SBarry Smith } 663b2f1ab58SBarry Smith 664b2f1ab58SBarry Smith if (do_values) { ierr = PetscFree(matrix_A->values);CHKERRQ(ierr);} 665b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->icols);CHKERRQ(ierr); 666b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->row_nnz);CHKERRQ(ierr); 667b2f1ab58SBarry Smith 6682205254eSKarl Rupp if (do_values) matrix_A->values = vals; 669b2f1ab58SBarry Smith matrix_A->icols = icols; 670b2f1ab58SBarry Smith matrix_A->row_nnz = row_nnz; 671b2f1ab58SBarry Smith PetscFunctionReturn(0); 672b2f1ab58SBarry Smith } 673b2f1ab58SBarry Smith 674b2f1ab58SBarry Smith /* 675b2f1ab58SBarry Smith spbas_apply_reordering_cols: 676b2f1ab58SBarry Smith apply the given reordering to the columns: matrix_A(:,perm) = matrix_A; 677b2f1ab58SBarry Smith */ 678b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_cols(spbas_matrix *matrix_A,const PetscInt *permutation) 679b2f1ab58SBarry Smith { 680b2f1ab58SBarry Smith PetscInt i,j; 681b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 682b2f1ab58SBarry Smith PetscInt row_nnz; 683b2f1ab58SBarry Smith PetscInt *icols; 684ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 6850298fd71SBarry Smith PetscScalar *vals = NULL; 686b2f1ab58SBarry Smith PetscErrorCode ierr; 687b2f1ab58SBarry Smith 688b2f1ab58SBarry Smith PetscFunctionBegin; 689e32f2f54SBarry 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"); 690b2f1ab58SBarry Smith 6916322f4bdSBarry Smith for (i=0; i<nrows; i++) { 692b2f1ab58SBarry Smith icols = matrix_A->icols[i]; 693b2f1ab58SBarry Smith row_nnz = matrix_A->row_nnz[i]; 6942205254eSKarl Rupp if (do_values) vals = matrix_A->values[i]; 695b2f1ab58SBarry Smith 6966322f4bdSBarry Smith for (j=0; j<row_nnz; j++) { 697b2f1ab58SBarry Smith icols[j] = permutation[i+icols[j]]-i; 698b2f1ab58SBarry Smith } 699b2f1ab58SBarry Smith ierr = spbas_mergesort(row_nnz, icols, vals);CHKERRQ(ierr); 700b2f1ab58SBarry Smith } 701b2f1ab58SBarry Smith PetscFunctionReturn(0); 702b2f1ab58SBarry Smith } 703b2f1ab58SBarry Smith 704b2f1ab58SBarry Smith /* 705b2f1ab58SBarry Smith spbas_apply_reordering: 706b2f1ab58SBarry Smith apply the given reordering: matrix_A(perm,perm) = matrix_A; 707b2f1ab58SBarry Smith */ 708b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering(spbas_matrix *matrix_A, const PetscInt *permutation, const PetscInt * inv_perm) 709b2f1ab58SBarry Smith { 710b2f1ab58SBarry Smith PetscErrorCode ierr; 7115fd66863SKarl Rupp 712b2f1ab58SBarry Smith PetscFunctionBegin; 713b2f1ab58SBarry Smith ierr = spbas_apply_reordering_rows(matrix_A, inv_perm);CHKERRQ(ierr); 714b2f1ab58SBarry Smith ierr = spbas_apply_reordering_cols(matrix_A, permutation);CHKERRQ(ierr); 715b2f1ab58SBarry Smith PetscFunctionReturn(0); 716b2f1ab58SBarry Smith } 717b2f1ab58SBarry Smith 718b2f1ab58SBarry Smith PetscErrorCode spbas_pattern_only(PetscInt nrows, PetscInt ncols, PetscInt *ai, PetscInt *aj, spbas_matrix * result) 719b2f1ab58SBarry Smith { 720b2f1ab58SBarry Smith spbas_matrix retval; 721b2f1ab58SBarry Smith PetscInt i, j, i0, r_nnz; 722b2f1ab58SBarry Smith PetscErrorCode ierr; 723b2f1ab58SBarry Smith 724b2f1ab58SBarry Smith PetscFunctionBegin; 725c328eeadSBarry Smith /* Copy input values */ 726b2f1ab58SBarry Smith retval.nrows = nrows; 727b2f1ab58SBarry Smith retval.ncols = ncols; 728b2f1ab58SBarry Smith retval.nnz = ai[nrows]; 729b2f1ab58SBarry Smith 730b2f1ab58SBarry Smith retval.block_data = PETSC_TRUE; 731b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 732b2f1ab58SBarry Smith 733c328eeadSBarry Smith /* Allocate output matrix */ 734b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(&retval, PETSC_FALSE);CHKERRQ(ierr); 7352205254eSKarl Rupp for (i=0; i<nrows; i++) retval.row_nnz[i] = ai[i+1]-ai[i]; 736b2f1ab58SBarry Smith ierr = spbas_allocate_data(&retval);CHKERRQ(ierr); 737c328eeadSBarry Smith /* Copy the structure */ 7386322f4bdSBarry Smith for (i = 0; i<retval.nrows; i++) { 739b2f1ab58SBarry Smith i0 = ai[i]; 740b2f1ab58SBarry Smith r_nnz = ai[i+1]-i0; 741b2f1ab58SBarry Smith 7426322f4bdSBarry Smith for (j=0; j<r_nnz; j++) { 743b2f1ab58SBarry Smith retval.icols[i][j] = aj[i0+j]-i; 744b2f1ab58SBarry Smith } 745b2f1ab58SBarry Smith } 746b2f1ab58SBarry Smith *result = retval; 747b2f1ab58SBarry Smith PetscFunctionReturn(0); 748b2f1ab58SBarry Smith } 749b2f1ab58SBarry Smith 750b2f1ab58SBarry Smith /* 751b2f1ab58SBarry Smith spbas_mark_row_power: 752b2f1ab58SBarry Smith Mark the columns in row 'row' which are nonzero in 753b2f1ab58SBarry Smith matrix^2log(marker). 754b2f1ab58SBarry Smith */ 755be332245SKarl Rupp PetscErrorCode spbas_mark_row_power(PetscInt *iwork, /* marker-vector */ 756c328eeadSBarry Smith PetscInt row, /* row for which the columns are marked */ 757c328eeadSBarry Smith spbas_matrix * in_matrix, /* matrix for which the power is being calculated */ 758c328eeadSBarry Smith PetscInt marker, /* marker-value: 2^power */ 759c328eeadSBarry Smith PetscInt minmrk, /* lower bound for marked points */ 760c328eeadSBarry Smith PetscInt maxmrk) /* upper bound for marked points */ 761b2f1ab58SBarry Smith { 762b2f1ab58SBarry Smith PetscErrorCode ierr; 763b2f1ab58SBarry Smith PetscInt i,j, nnz; 764b2f1ab58SBarry Smith 765b2f1ab58SBarry Smith PetscFunctionBegin; 766b2f1ab58SBarry Smith nnz = in_matrix->row_nnz[row]; 767b2f1ab58SBarry Smith 768c328eeadSBarry Smith /* For higher powers, call this function recursively */ 7696322f4bdSBarry Smith if (marker>1) { 7706322f4bdSBarry Smith for (i=0; i<nnz; i++) { 771b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 7726322f4bdSBarry Smith if (minmrk<=j && j<maxmrk && iwork[j] < marker) { 77371d55d18SBarry Smith ierr = spbas_mark_row_power(iwork, row + in_matrix->icols[row][i],in_matrix, marker/2,minmrk,maxmrk);CHKERRQ(ierr); 774b2f1ab58SBarry Smith iwork[j] |= marker; 775b2f1ab58SBarry Smith } 776b2f1ab58SBarry Smith } 7776322f4bdSBarry Smith } else { 778c328eeadSBarry Smith /* Mark the columns reached */ 7796322f4bdSBarry Smith for (i=0; i<nnz; i++) { 780b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 7812205254eSKarl Rupp if (minmrk<=j && j<maxmrk) iwork[j] |= 1; 782b2f1ab58SBarry Smith } 783b2f1ab58SBarry Smith } 784b2f1ab58SBarry Smith PetscFunctionReturn(0); 785b2f1ab58SBarry Smith } 786b2f1ab58SBarry Smith 787b2f1ab58SBarry Smith /* 788b2f1ab58SBarry Smith spbas_power 789b2f1ab58SBarry Smith Calculate sparseness patterns for incomplete Cholesky decompositions 790b2f1ab58SBarry Smith of a given order: (almost) all nonzeros of the matrix^(order+1) which 791b2f1ab58SBarry Smith are inside the band width are found and stored in the output sparseness 792b2f1ab58SBarry Smith pattern. 793b2f1ab58SBarry Smith */ 794b2f1ab58SBarry Smith PetscErrorCode spbas_power(spbas_matrix in_matrix,PetscInt power, spbas_matrix * result) 795b2f1ab58SBarry Smith { 796b2f1ab58SBarry Smith spbas_matrix retval; 797b2f1ab58SBarry Smith PetscInt nrows = in_matrix.nrows; 798b2f1ab58SBarry Smith PetscInt ncols = in_matrix.ncols; 799b2f1ab58SBarry Smith PetscInt i, j, kend; 800b2f1ab58SBarry Smith PetscInt nnz, inz; 801b2f1ab58SBarry Smith PetscInt *iwork; 802b2f1ab58SBarry Smith PetscInt marker; 803b2f1ab58SBarry Smith PetscInt maxmrk=0; 804b2f1ab58SBarry Smith PetscErrorCode ierr; 805b2f1ab58SBarry Smith 806b2f1ab58SBarry Smith PetscFunctionBegin; 807e32f2f54SBarry 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"); 808e32f2f54SBarry Smith if (ncols != nrows) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Dimension error\n"); 809e32f2f54SBarry Smith if (in_matrix.values) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input array must be sparseness pattern (no values)"); 810e32f2f54SBarry Smith if (power<=0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Power must be 1 or up"); 811b2f1ab58SBarry Smith 812c328eeadSBarry Smith /* Copy input values*/ 813b2f1ab58SBarry Smith retval.nrows = ncols; 814b2f1ab58SBarry Smith retval.ncols = nrows; 815b2f1ab58SBarry Smith retval.nnz = 0; 816b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 817b2f1ab58SBarry Smith retval.block_data = PETSC_FALSE; 818b2f1ab58SBarry Smith 819c328eeadSBarry Smith /* Allocate sparseness pattern */ 82071d55d18SBarry Smith ierr = spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 821b2f1ab58SBarry Smith 822580bdb30SBarry Smith /* Allocate marker array: note sure the max needed so use the max of the two */ 823580bdb30SBarry Smith ierr = PetscCalloc1(PetscMax(ncols,nrows), &iwork);CHKERRQ(ierr); 824b2f1ab58SBarry Smith 825c328eeadSBarry Smith /* Calculate marker values */ 8262205254eSKarl Rupp marker = 1; for (i=1; i<power; i++) marker*=2; 827b2f1ab58SBarry Smith 8286322f4bdSBarry Smith for (i=0; i<nrows; i++) { 829c328eeadSBarry Smith /* Calculate the pattern for each row */ 830b2f1ab58SBarry Smith 831b2f1ab58SBarry Smith nnz = in_matrix.row_nnz[i]; 832b2f1ab58SBarry Smith kend = i+in_matrix.icols[i][nnz-1]; 8332205254eSKarl Rupp if (maxmrk<=kend) maxmrk=kend+1; 83471d55d18SBarry Smith ierr = spbas_mark_row_power(iwork, i, &in_matrix, marker, i, maxmrk);CHKERRQ(ierr); 835b2f1ab58SBarry Smith 836c328eeadSBarry Smith /* Count the columns*/ 837b2f1ab58SBarry Smith nnz = 0; 8382205254eSKarl Rupp for (j=i; j<maxmrk; j++) nnz+= (iwork[j]!=0); 839b2f1ab58SBarry Smith 840c328eeadSBarry Smith /* Allocate the column indices */ 841b2f1ab58SBarry Smith retval.row_nnz[i] = nnz; 842785e854fSJed Brown ierr = PetscMalloc1(nnz,&retval.icols[i]);CHKERRQ(ierr); 843b2f1ab58SBarry Smith 844c328eeadSBarry Smith /* Administrate the column indices */ 845b2f1ab58SBarry Smith inz = 0; 8466322f4bdSBarry Smith for (j=i; j<maxmrk; j++) { 8476322f4bdSBarry Smith if (iwork[j]) { 848b2f1ab58SBarry Smith retval.icols[i][inz] = j-i; 849b2f1ab58SBarry Smith inz++; 850b2f1ab58SBarry Smith iwork[j] = 0; 851b2f1ab58SBarry Smith } 852b2f1ab58SBarry Smith } 853b2f1ab58SBarry Smith retval.nnz += nnz; 854b2f1ab58SBarry Smith }; 855b2f1ab58SBarry Smith ierr = PetscFree(iwork);CHKERRQ(ierr); 856b2f1ab58SBarry Smith *result = retval; 857b2f1ab58SBarry Smith PetscFunctionReturn(0); 858b2f1ab58SBarry Smith } 859b2f1ab58SBarry Smith 860b2f1ab58SBarry Smith /* 861b2f1ab58SBarry Smith spbas_keep_upper: 862b2f1ab58SBarry Smith remove the lower part of the matrix: keep the upper part 863b2f1ab58SBarry Smith */ 864b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper(spbas_matrix * inout_matrix) 865b2f1ab58SBarry Smith { 866b2f1ab58SBarry Smith PetscInt i, j; 867b2f1ab58SBarry Smith PetscInt jstart; 868b2f1ab58SBarry Smith 869b2f1ab58SBarry Smith PetscFunctionBegin; 870e32f2f54SBarry Smith if (inout_matrix->block_data) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Not yet for block data matrices\n"); 8716322f4bdSBarry Smith for (i=0; i<inout_matrix->nrows; i++) { 8726322f4bdSBarry Smith for (jstart=0; (jstart<inout_matrix->row_nnz[i]) && (inout_matrix->icols[i][jstart]<0); jstart++) {} 8736322f4bdSBarry Smith if (jstart>0) { 8746322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 8756322f4bdSBarry Smith inout_matrix->icols[i][j] = inout_matrix->icols[i][j+jstart]; 876b2f1ab58SBarry Smith } 877b2f1ab58SBarry Smith 8786322f4bdSBarry Smith if (inout_matrix->values) { 8796322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 8806322f4bdSBarry Smith inout_matrix->values[i][j] = inout_matrix->values[i][j+jstart]; 881b2f1ab58SBarry Smith } 882b2f1ab58SBarry Smith } 883b2f1ab58SBarry Smith 884b2f1ab58SBarry Smith inout_matrix->row_nnz[i] -= jstart; 885b2f1ab58SBarry Smith 8866322f4bdSBarry Smith inout_matrix->icols[i] = (PetscInt*) realloc((void*) inout_matrix->icols[i], inout_matrix->row_nnz[i]*sizeof(PetscInt)); 887b2f1ab58SBarry Smith 8886322f4bdSBarry Smith if (inout_matrix->values) { 8896322f4bdSBarry Smith inout_matrix->values[i] = (PetscScalar*) realloc((void*) inout_matrix->values[i], inout_matrix->row_nnz[i]*sizeof(PetscScalar)); 890b2f1ab58SBarry Smith } 891b2f1ab58SBarry Smith inout_matrix->nnz -= jstart; 892b2f1ab58SBarry Smith } 893b2f1ab58SBarry Smith } 894b2f1ab58SBarry Smith PetscFunctionReturn(0); 895b2f1ab58SBarry Smith } 896