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 PetscInt nrows = result->nrows; 64b2f1ab58SBarry Smith PetscInt col_idx_type = result->col_idx_type; 65b2f1ab58SBarry Smith 66b2f1ab58SBarry Smith PetscFunctionBegin; 67c328eeadSBarry Smith /* Allocate sparseness pattern */ 689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows,&result->row_nnz)); 699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows,&result->icols)); 70b2f1ab58SBarry Smith 71c328eeadSBarry Smith /* If offsets are given wrt an array, create array */ 724e1ff37aSBarry Smith if (col_idx_type == SPBAS_OFFSET_ARRAY) { 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows,&result->icol0)); 744e1ff37aSBarry Smith } else { 750298fd71SBarry Smith result->icol0 = NULL; 76b2f1ab58SBarry Smith } 77b2f1ab58SBarry Smith 78c328eeadSBarry Smith /* If values are given, allocate values array */ 794e1ff37aSBarry Smith if (do_values) { 809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows,&result->values)); 814e1ff37aSBarry Smith } else { 820298fd71SBarry Smith result->values = NULL; 83b2f1ab58SBarry Smith } 84b2f1ab58SBarry Smith PetscFunctionReturn(0); 85b2f1ab58SBarry Smith } 86b2f1ab58SBarry Smith 87b2f1ab58SBarry Smith /* 88b2f1ab58SBarry Smith spbas_allocate_data: 89b2f1ab58SBarry Smith in case of block_data: 90b2f1ab58SBarry Smith Allocate the data arrays alloc_icol and optionally alloc_val, 91b2f1ab58SBarry Smith set appropriate pointers from icols and values; 92b2f1ab58SBarry Smith in case of !block_data: 93b2f1ab58SBarry Smith Allocate the arrays icols[i] and optionally values[i] 94b2f1ab58SBarry Smith */ 95b2f1ab58SBarry Smith PetscErrorCode spbas_allocate_data(spbas_matrix * result) 96b2f1ab58SBarry Smith { 97b2f1ab58SBarry Smith PetscInt i; 98b2f1ab58SBarry Smith PetscInt nnz = result->nnz; 99b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 100b2f1ab58SBarry Smith PetscInt r_nnz; 1016c4ed002SBarry Smith PetscBool do_values = (result->values) ? PETSC_TRUE : PETSC_FALSE; 102ace3abfcSBarry Smith PetscBool block_data = result->block_data; 103b2f1ab58SBarry Smith 104b2f1ab58SBarry Smith PetscFunctionBegin; 1054e1ff37aSBarry Smith if (block_data) { 106c328eeadSBarry Smith /* Allocate the column number array and point to it */ 107b2f1ab58SBarry Smith result->n_alloc_icol = nnz; 1082205254eSKarl Rupp 1099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz, &result->alloc_icol)); 1102205254eSKarl Rupp 111b2f1ab58SBarry Smith result->icols[0] = result->alloc_icol; 1124e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 113b2f1ab58SBarry Smith result->icols[i] = result->icols[i-1] + result->row_nnz[i-1]; 114b2f1ab58SBarry Smith } 115b2f1ab58SBarry Smith 116c328eeadSBarry Smith /* Allocate the value array and point to it */ 1174e1ff37aSBarry Smith if (do_values) { 118b2f1ab58SBarry Smith result->n_alloc_val = nnz; 1192205254eSKarl Rupp 1209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz, &result->alloc_val)); 1212205254eSKarl Rupp 122b2f1ab58SBarry Smith result->values[0] = result->alloc_val; 1234e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 124b2f1ab58SBarry Smith result->values[i] = result->values[i-1] + result->row_nnz[i-1]; 125b2f1ab58SBarry Smith } 126b2f1ab58SBarry Smith } 1274e1ff37aSBarry Smith } else { 1284e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 129b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 1309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(r_nnz, &result->icols[i])); 131b2f1ab58SBarry Smith } 1324e1ff37aSBarry Smith if (do_values) { 1334e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 134b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 1359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(r_nnz, &result->values[i])); 136b2f1ab58SBarry Smith } 137b2f1ab58SBarry Smith } 138b2f1ab58SBarry Smith } 139b2f1ab58SBarry Smith PetscFunctionReturn(0); 140b2f1ab58SBarry Smith } 141b2f1ab58SBarry Smith 142b2f1ab58SBarry Smith /* 143b2f1ab58SBarry Smith spbas_row_order_icol 144b2f1ab58SBarry Smith determine if row i1 should come 145b2f1ab58SBarry Smith + before row i2 in the sorted rows (return -1), 146b2f1ab58SBarry Smith + after (return 1) 147b2f1ab58SBarry Smith + is identical (return 0). 148b2f1ab58SBarry Smith */ 149b2f1ab58SBarry Smith int spbas_row_order_icol(PetscInt i1, PetscInt i2, PetscInt *irow_in, PetscInt *icol_in,PetscInt col_idx_type) 150b2f1ab58SBarry Smith { 151b2f1ab58SBarry Smith PetscInt j; 152b2f1ab58SBarry Smith PetscInt nnz1 = irow_in[i1+1] - irow_in[i1]; 153b2f1ab58SBarry Smith PetscInt nnz2 = irow_in[i2+1] - irow_in[i2]; 154b2f1ab58SBarry Smith PetscInt * icol1 = &icol_in[irow_in[i1]]; 155b2f1ab58SBarry Smith PetscInt * icol2 = &icol_in[irow_in[i2]]; 156b2f1ab58SBarry Smith 1572205254eSKarl Rupp if (nnz1<nnz2) return -1; 1582205254eSKarl Rupp if (nnz1>nnz2) return 1; 159b2f1ab58SBarry Smith 1604e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 1614e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 1622205254eSKarl Rupp if (icol1[j]< icol2[j]) return -1; 1632205254eSKarl Rupp if (icol1[j]> icol2[j]) return 1; 164b2f1ab58SBarry Smith } 1654e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 1664e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 1672205254eSKarl Rupp if (icol1[j]-i1< icol2[j]-i2) return -1; 1682205254eSKarl Rupp if (icol1[j]-i1> icol2[j]-i2) return 1; 169b2f1ab58SBarry Smith } 1704e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 1714e1ff37aSBarry Smith for (j=1; j<nnz1; j++) { 1722205254eSKarl Rupp if (icol1[j]-icol1[0] < icol2[j]-icol2[0]) return -1; 1732205254eSKarl Rupp if (icol1[j]-icol1[0] > icol2[j]-icol2[0]) return 1; 174b2f1ab58SBarry Smith } 175b2f1ab58SBarry Smith } 176b2f1ab58SBarry Smith return 0; 177b2f1ab58SBarry Smith } 178b2f1ab58SBarry Smith 179b2f1ab58SBarry Smith /* 180b2f1ab58SBarry Smith spbas_mergesort_icols: 181b2f1ab58SBarry Smith return a sorting of the rows in which identical sparseness patterns are 182b2f1ab58SBarry Smith next to each other 183b2f1ab58SBarry Smith */ 184b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort_icols(PetscInt nrows, PetscInt * irow_in, PetscInt * icol_in,PetscInt col_idx_type, PetscInt *isort) 185b2f1ab58SBarry Smith { 186c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 187c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 188c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 189c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 190c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 191c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 192c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 193b2f1ab58SBarry Smith 194b2f1ab58SBarry Smith PetscFunctionBegin; 1959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows,&ialloc)); 196b2f1ab58SBarry Smith 197b2f1ab58SBarry Smith ihlp1 = ialloc; 198b2f1ab58SBarry Smith ihlp2 = isort; 199b2f1ab58SBarry Smith 200c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 2014e1ff37aSBarry Smith for (istep=1; istep<nrows; istep*=2) { 202c328eeadSBarry Smith /* 203c328eeadSBarry Smith Combine sorted parts 204c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 205c328eeadSBarry Smith of ihlp2 and vhlp2 206c328eeadSBarry Smith 207c328eeadSBarry Smith into one sorted part 208c328eeadSBarry Smith istart:istart+2*istep-1 209c328eeadSBarry Smith of ihlp1 and vhlp1 210c328eeadSBarry Smith */ 2114e1ff37aSBarry Smith for (istart=0; istart<nrows; istart+=2*istep) { 212c328eeadSBarry Smith /* Set counters and bound array part endings */ 2132205254eSKarl Rupp i1=istart; i1end = i1+istep; if (i1end>nrows) i1end=nrows; 2142205254eSKarl Rupp i2=istart+istep; i2end = i2+istep; if (i2end>nrows) i2end=nrows; 215b2f1ab58SBarry Smith 216c328eeadSBarry Smith /* Merge the two array parts */ 2174e1ff37aSBarry Smith for (i=istart; i<i2end; i++) { 2184e1ff37aSBarry Smith if (i1<i1end && i2<i2end && spbas_row_order_icol(ihlp2[i1], ihlp2[i2], irow_in, icol_in, col_idx_type) < 0) { 219b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 220b2f1ab58SBarry Smith i1++; 2214e1ff37aSBarry Smith } else if (i2<i2end) { 222b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 223b2f1ab58SBarry Smith i2++; 2244e1ff37aSBarry Smith } else { 225b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 226b2f1ab58SBarry Smith i1++; 227b2f1ab58SBarry Smith } 228b2f1ab58SBarry Smith } 229b2f1ab58SBarry Smith } 230b2f1ab58SBarry Smith 231c328eeadSBarry Smith /* Swap the two array sets */ 232b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 233b2f1ab58SBarry Smith } 234b2f1ab58SBarry Smith 235c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 2364e1ff37aSBarry Smith if (ihlp2 != isort) { 2372205254eSKarl Rupp for (i=0; i<nrows; i++) isort[i] = ihlp2[i]; 238b2f1ab58SBarry Smith } 2399566063dSJacob Faibussowitsch PetscCall(PetscFree(ialloc)); 240b2f1ab58SBarry Smith PetscFunctionReturn(0); 241b2f1ab58SBarry Smith } 242b2f1ab58SBarry Smith 243b2f1ab58SBarry Smith /* 244b2f1ab58SBarry Smith spbas_compress_pattern: 245b2f1ab58SBarry Smith calculate a compressed sparseness pattern for a sparseness pattern 246b2f1ab58SBarry Smith given in compressed row storage. The compressed sparseness pattern may 247b2f1ab58SBarry Smith require (much) less memory. 248b2f1ab58SBarry Smith */ 2494e1ff37aSBarry 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) 250b2f1ab58SBarry Smith { 251b2f1ab58SBarry Smith PetscInt nnz = irow_in[nrows]; 2523dfa2556SBarry Smith size_t mem_orig = (nrows + nnz) * sizeof(PetscInt); 2533dfa2556SBarry Smith size_t mem_compressed; 254b2f1ab58SBarry Smith PetscInt *isort; 255b2f1ab58SBarry Smith PetscInt *icols; 256b2f1ab58SBarry Smith PetscInt row_nnz; 257b2f1ab58SBarry Smith PetscInt *ipoint; 258ace3abfcSBarry Smith PetscBool *used; 259b2f1ab58SBarry Smith PetscInt ptr; 260b2f1ab58SBarry Smith PetscInt i,j; 261ace3abfcSBarry Smith const PetscBool no_values = PETSC_FALSE; 262b2f1ab58SBarry Smith 263b2f1ab58SBarry Smith PetscFunctionBegin; 264c328eeadSBarry Smith /* Allocate the structure of the new matrix */ 265b2f1ab58SBarry Smith B->nrows = nrows; 266b2f1ab58SBarry Smith B->ncols = ncols; 267b2f1ab58SBarry Smith B->nnz = nnz; 268b2f1ab58SBarry Smith B->col_idx_type = col_idx_type; 269b2f1ab58SBarry Smith B->block_data = PETSC_TRUE; 2702205254eSKarl Rupp 2719566063dSJacob Faibussowitsch PetscCall(spbas_allocate_pattern(B, no_values)); 272b2f1ab58SBarry Smith 273c328eeadSBarry Smith /* When using an offset array, set it */ 2744e1ff37aSBarry Smith if (col_idx_type==SPBAS_OFFSET_ARRAY) { 2752205254eSKarl Rupp for (i=0; i<nrows; i++) B->icol0[i] = icol_in[irow_in[i]]; 276b2f1ab58SBarry Smith } 277b2f1ab58SBarry Smith 278c328eeadSBarry Smith /* Allocate the ordering for the rows */ 2799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows,&isort)); 2809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows,&ipoint)); 2819566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nrows,&used)); 282b2f1ab58SBarry Smith 2834e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 284b2f1ab58SBarry Smith B->row_nnz[i] = irow_in[i+1]-irow_in[i]; 285b2f1ab58SBarry Smith isort[i] = i; 286b2f1ab58SBarry Smith ipoint[i] = i; 287b2f1ab58SBarry Smith } 288b2f1ab58SBarry Smith 289c328eeadSBarry Smith /* Sort the rows so that identical columns will be next to each other */ 2909566063dSJacob Faibussowitsch PetscCall(spbas_mergesort_icols(nrows, irow_in, icol_in, col_idx_type, isort)); 2919566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"Rows have been sorted for patterns\n")); 292b2f1ab58SBarry Smith 293c328eeadSBarry Smith /* Replace identical rows with the first one in the list */ 2944e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 2954e1ff37aSBarry Smith if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) { 296b2f1ab58SBarry Smith ipoint[isort[i]] = ipoint[isort[i-1]]; 297b2f1ab58SBarry Smith } 298b2f1ab58SBarry Smith } 299b2f1ab58SBarry Smith 300c328eeadSBarry Smith /* Collect the rows which are used*/ 3012205254eSKarl Rupp for (i=0; i<nrows; i++) used[ipoint[i]] = PETSC_TRUE; 302b2f1ab58SBarry Smith 303c328eeadSBarry Smith /* Calculate needed memory */ 304b2f1ab58SBarry Smith B->n_alloc_icol = 0; 3054e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 3062205254eSKarl Rupp if (used[i]) B->n_alloc_icol += B->row_nnz[i]; 307b2f1ab58SBarry Smith } 3089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(B->n_alloc_icol,&B->alloc_icol)); 309b2f1ab58SBarry Smith 310c328eeadSBarry Smith /* Fill in the diagonal offsets for the rows which store their own data */ 311b2f1ab58SBarry Smith ptr = 0; 3124e1ff37aSBarry Smith for (i=0; i<B->nrows; i++) { 3134e1ff37aSBarry Smith if (used[i]) { 314b2f1ab58SBarry Smith B->icols[i] = &B->alloc_icol[ptr]; 315b2f1ab58SBarry Smith icols = &icol_in[irow_in[i]]; 316b2f1ab58SBarry Smith row_nnz = B->row_nnz[i]; 3174e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 3184e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 319b2f1ab58SBarry Smith B->icols[i][j] = icols[j]; 320b2f1ab58SBarry Smith } 3214e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 3224e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 323b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-i; 324b2f1ab58SBarry Smith } 3254e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 3264e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 327b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-icols[0]; 328b2f1ab58SBarry Smith } 329b2f1ab58SBarry Smith } 330b2f1ab58SBarry Smith ptr += B->row_nnz[i]; 331b2f1ab58SBarry Smith } 332b2f1ab58SBarry Smith } 333b2f1ab58SBarry Smith 334c328eeadSBarry Smith /* Point to the right places for all data */ 3354e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 336b2f1ab58SBarry Smith B->icols[i] = B->icols[ipoint[i]]; 337b2f1ab58SBarry Smith } 3389566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"Row patterns have been compressed\n")); 3399566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL," (%g nonzeros per row)\n", (double) ((PetscReal) nnz / (PetscReal) nrows))); 340b2f1ab58SBarry Smith 3419566063dSJacob Faibussowitsch PetscCall(PetscFree(isort)); 3429566063dSJacob Faibussowitsch PetscCall(PetscFree(used)); 3439566063dSJacob Faibussowitsch PetscCall(PetscFree(ipoint)); 344b2f1ab58SBarry Smith 345b2f1ab58SBarry Smith mem_compressed = spbas_memory_requirement(*B); 3464e1ff37aSBarry Smith *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig; 347b2f1ab58SBarry Smith PetscFunctionReturn(0); 348b2f1ab58SBarry Smith } 349b2f1ab58SBarry Smith 350b2f1ab58SBarry Smith /* 351b2f1ab58SBarry Smith spbas_incomplete_cholesky 352b2f1ab58SBarry Smith Incomplete Cholesky decomposition 353b2f1ab58SBarry Smith */ 354c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas_cholesky.h> 355b2f1ab58SBarry Smith 356b2f1ab58SBarry Smith /* 357b2f1ab58SBarry Smith spbas_delete : de-allocate the arrays owned by this matrix 358b2f1ab58SBarry Smith */ 359b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix matrix) 360b2f1ab58SBarry Smith { 361b2f1ab58SBarry Smith PetscInt i; 3629ccc4a9bSBarry Smith 363b2f1ab58SBarry Smith PetscFunctionBegin; 3649ccc4a9bSBarry Smith if (matrix.block_data) { 3659566063dSJacob Faibussowitsch PetscCall(PetscFree(matrix.alloc_icol)); 3669566063dSJacob Faibussowitsch if (matrix.values) PetscCall(PetscFree(matrix.alloc_val)); 3679ccc4a9bSBarry Smith } else { 3689566063dSJacob Faibussowitsch for (i=0; i<matrix.nrows; i++) PetscCall(PetscFree(matrix.icols[i])); 3699566063dSJacob Faibussowitsch PetscCall(PetscFree(matrix.icols)); 3709ccc4a9bSBarry Smith if (matrix.values) { 3719566063dSJacob Faibussowitsch for (i=0; i<matrix.nrows; i++) PetscCall(PetscFree(matrix.values[i])); 372b2f1ab58SBarry Smith } 373b2f1ab58SBarry Smith } 374b2f1ab58SBarry Smith 3759566063dSJacob Faibussowitsch PetscCall(PetscFree(matrix.row_nnz)); 3769566063dSJacob Faibussowitsch PetscCall(PetscFree(matrix.icols)); 3779566063dSJacob Faibussowitsch if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) PetscCall(PetscFree(matrix.icol0)); 3789566063dSJacob Faibussowitsch PetscCall(PetscFree(matrix.values)); 379b2f1ab58SBarry Smith PetscFunctionReturn(0); 380b2f1ab58SBarry Smith } 381b2f1ab58SBarry Smith 382b2f1ab58SBarry Smith /* 383b2f1ab58SBarry Smith spbas_matrix_to_crs: 384b2f1ab58SBarry Smith Convert an spbas_matrix to compessed row storage 385b2f1ab58SBarry Smith */ 386b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out) 387b2f1ab58SBarry Smith { 388b2f1ab58SBarry Smith PetscInt nrows = matrix_A.nrows; 389b2f1ab58SBarry Smith PetscInt nnz = matrix_A.nnz; 390b2f1ab58SBarry Smith PetscInt i,j,r_nnz,i0; 391b2f1ab58SBarry Smith PetscInt *irow; 392b2f1ab58SBarry Smith PetscInt *icol; 393b2f1ab58SBarry Smith PetscInt *icol_A; 394b2f1ab58SBarry Smith MatScalar *val; 395b2f1ab58SBarry Smith PetscScalar *val_A; 396b2f1ab58SBarry Smith PetscInt col_idx_type = matrix_A.col_idx_type; 397ace3abfcSBarry Smith PetscBool do_values = matrix_A.values ? PETSC_TRUE : PETSC_FALSE; 398b2f1ab58SBarry Smith 399b2f1ab58SBarry Smith PetscFunctionBegin; 4009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows+1, &irow)); 4019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz, &icol)); 4029ccc4a9bSBarry Smith *icol_out = icol; 4039ccc4a9bSBarry Smith *irow_out = irow; 4049ccc4a9bSBarry Smith if (do_values) { 4059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz, &val)); 406b2f1ab58SBarry Smith *val_out = val; *icol_out = icol; *irow_out=irow; 407b2f1ab58SBarry Smith } 408b2f1ab58SBarry Smith 409b2f1ab58SBarry Smith irow[0]=0; 4109ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 411b2f1ab58SBarry Smith r_nnz = matrix_A.row_nnz[i]; 412b2f1ab58SBarry Smith i0 = irow[i]; 413b2f1ab58SBarry Smith irow[i+1] = i0 + r_nnz; 414b2f1ab58SBarry Smith icol_A = matrix_A.icols[i]; 415b2f1ab58SBarry Smith 4169ccc4a9bSBarry Smith if (do_values) { 417b2f1ab58SBarry Smith val_A = matrix_A.values[i]; 4189ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 419b2f1ab58SBarry Smith icol[i0+j] = icol_A[j]; 420b2f1ab58SBarry Smith val[i0+j] = val_A[j]; 421b2f1ab58SBarry Smith } 4229ccc4a9bSBarry Smith } else { 4232205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] = icol_A[j]; 424b2f1ab58SBarry Smith } 425b2f1ab58SBarry Smith 4269ccc4a9bSBarry Smith if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 4272205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] += i; 4282205254eSKarl Rupp } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 429b2f1ab58SBarry Smith i0 = matrix_A.icol0[i]; 4302205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] += i0; 431b2f1ab58SBarry Smith } 432b2f1ab58SBarry Smith } 433b2f1ab58SBarry Smith PetscFunctionReturn(0); 434b2f1ab58SBarry Smith } 435b2f1ab58SBarry Smith 436b2f1ab58SBarry Smith /* 437b2f1ab58SBarry Smith spbas_transpose 438b2f1ab58SBarry Smith return the transpose of a matrix 439b2f1ab58SBarry Smith */ 440b2f1ab58SBarry Smith PetscErrorCode spbas_transpose(spbas_matrix in_matrix, spbas_matrix * result) 441b2f1ab58SBarry Smith { 442b2f1ab58SBarry Smith PetscInt col_idx_type = in_matrix.col_idx_type; 443b2f1ab58SBarry Smith PetscInt nnz = in_matrix.nnz; 444b2f1ab58SBarry Smith PetscInt ncols = in_matrix.nrows; 445b2f1ab58SBarry Smith PetscInt nrows = in_matrix.ncols; 446b2f1ab58SBarry Smith PetscInt i,j,k; 447b2f1ab58SBarry Smith PetscInt r_nnz; 448b2f1ab58SBarry Smith PetscInt *irow; 4494efc9174SBarry Smith PetscInt icol0 = 0; 450b2f1ab58SBarry Smith PetscScalar * val; 4514e1ff37aSBarry Smith 452b2f1ab58SBarry Smith PetscFunctionBegin; 453c328eeadSBarry Smith /* Copy input values */ 454b2f1ab58SBarry Smith result->nrows = nrows; 455b2f1ab58SBarry Smith result->ncols = ncols; 456b2f1ab58SBarry Smith result->nnz = nnz; 457b2f1ab58SBarry Smith result->col_idx_type = SPBAS_COLUMN_NUMBERS; 458b2f1ab58SBarry Smith result->block_data = PETSC_TRUE; 459b2f1ab58SBarry Smith 460c328eeadSBarry Smith /* Allocate sparseness pattern */ 4619566063dSJacob Faibussowitsch PetscCall(spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE)); 462b2f1ab58SBarry Smith 463c328eeadSBarry Smith /* Count the number of nonzeros in each row */ 4642205254eSKarl Rupp for (i = 0; i<nrows; i++) result->row_nnz[i] = 0; 465b2f1ab58SBarry Smith 4669ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 467b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 468b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 4699ccc4a9bSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 4702205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[irow[j]]++; 4719ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 4722205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[i+irow[j]]++; 4739ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 474b2f1ab58SBarry Smith icol0=in_matrix.icol0[i]; 4752205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[icol0+irow[j]]++; 476b2f1ab58SBarry Smith } 477b2f1ab58SBarry Smith } 478b2f1ab58SBarry Smith 479c328eeadSBarry Smith /* Set the pointers to the data */ 4809566063dSJacob Faibussowitsch PetscCall(spbas_allocate_data(result)); 481b2f1ab58SBarry Smith 482c328eeadSBarry Smith /* Reset the number of nonzeros in each row */ 4832205254eSKarl Rupp for (i = 0; i<nrows; i++) result->row_nnz[i] = 0; 484b2f1ab58SBarry Smith 485c328eeadSBarry Smith /* Fill the data arrays */ 4869ccc4a9bSBarry Smith if (in_matrix.values) { 4879ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 488b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 489b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 490b2f1ab58SBarry Smith val = in_matrix.values[i]; 491b2f1ab58SBarry Smith 4922205254eSKarl Rupp if (col_idx_type == SPBAS_COLUMN_NUMBERS) icol0 = 0; 4932205254eSKarl Rupp else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0 = i; 4942205254eSKarl Rupp else if (col_idx_type == SPBAS_OFFSET_ARRAY) icol0 = in_matrix.icol0[i]; 4959ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 496b2f1ab58SBarry Smith k = icol0 + irow[j]; 497b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 498b2f1ab58SBarry Smith result->values[k][result->row_nnz[k]] = val[j]; 499b2f1ab58SBarry Smith result->row_nnz[k]++; 500b2f1ab58SBarry Smith } 501b2f1ab58SBarry Smith } 5029ccc4a9bSBarry Smith } else { 5039ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 504b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 505b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 506b2f1ab58SBarry Smith 5072205254eSKarl Rupp if (col_idx_type == SPBAS_COLUMN_NUMBERS) icol0=0; 5082205254eSKarl Rupp else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0=i; 5092205254eSKarl Rupp else if (col_idx_type == SPBAS_OFFSET_ARRAY) icol0=in_matrix.icol0[i]; 510b2f1ab58SBarry Smith 5119ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 512b2f1ab58SBarry Smith k = icol0 + irow[j]; 513b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 514b2f1ab58SBarry Smith result->row_nnz[k]++; 515b2f1ab58SBarry Smith } 516b2f1ab58SBarry Smith } 517b2f1ab58SBarry Smith } 518b2f1ab58SBarry Smith PetscFunctionReturn(0); 519b2f1ab58SBarry Smith } 520b2f1ab58SBarry Smith 521b2f1ab58SBarry Smith /* 522b2f1ab58SBarry Smith spbas_mergesort 523b2f1ab58SBarry Smith 524bb42e86aSJed Brown mergesort for an array of integers and an array of associated 525b2f1ab58SBarry Smith reals 526b2f1ab58SBarry Smith 527b2f1ab58SBarry Smith on output, icol[0..nnz-1] is increasing; 528b2f1ab58SBarry Smith val[0..nnz-1] has undergone the same permutation as icol 529b2f1ab58SBarry Smith 5300298fd71SBarry Smith NB: val may be NULL: in that case, only the integers are sorted 531b2f1ab58SBarry Smith 532b2f1ab58SBarry Smith */ 533b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val) 534b2f1ab58SBarry Smith { 535c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 536c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 537c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 538c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 5390298fd71SBarry Smith PetscScalar *valloc=NULL; 540c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 541b2f1ab58SBarry Smith PetscScalar *vswap; 542c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 5430298fd71SBarry Smith PetscScalar *vhlp1=NULL; /* (arrays under construction) */ 544c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 5450298fd71SBarry Smith PetscScalar *vhlp2=NULL; 546b2f1ab58SBarry Smith 547362febeeSStefano Zampini PetscFunctionBegin; 5489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz,&ialloc)); 549b2f1ab58SBarry Smith ihlp1 = ialloc; 550b2f1ab58SBarry Smith ihlp2 = icol; 551b2f1ab58SBarry Smith 5529ccc4a9bSBarry Smith if (val) { 5539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz,&valloc)); 554b2f1ab58SBarry Smith vhlp1 = valloc; 555b2f1ab58SBarry Smith vhlp2 = val; 556b2f1ab58SBarry Smith } 557b2f1ab58SBarry Smith 558c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 5599ccc4a9bSBarry Smith for (istep=1; istep<nnz; istep*=2) { 560c328eeadSBarry Smith /* 561c328eeadSBarry Smith Combine sorted parts 562c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 563c328eeadSBarry Smith of ihlp2 and vhlp2 564c328eeadSBarry Smith 565c328eeadSBarry Smith into one sorted part 566c328eeadSBarry Smith istart:istart+2*istep-1 567c328eeadSBarry Smith of ihlp1 and vhlp1 568c328eeadSBarry Smith */ 5699ccc4a9bSBarry Smith for (istart=0; istart<nnz; istart+=2*istep) { 570c328eeadSBarry Smith /* Set counters and bound array part endings */ 5712205254eSKarl Rupp i1=istart; i1end = i1+istep; if (i1end>nnz) i1end=nnz; 5722205254eSKarl Rupp i2=istart+istep; i2end = i2+istep; if (i2end>nnz) i2end=nnz; 573b2f1ab58SBarry Smith 574c328eeadSBarry Smith /* Merge the two array parts */ 5759ccc4a9bSBarry Smith if (val) { 5769ccc4a9bSBarry Smith for (i=istart; i<i2end; i++) { 5779ccc4a9bSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 578b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 579b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 580b2f1ab58SBarry Smith i1++; 5819ccc4a9bSBarry Smith } else if (i2<i2end) { 582b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 583b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i2]; 584b2f1ab58SBarry Smith i2++; 5859ccc4a9bSBarry Smith } else { 586b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 587b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 588b2f1ab58SBarry Smith i1++; 589b2f1ab58SBarry Smith } 590b2f1ab58SBarry Smith } 5919ccc4a9bSBarry Smith } else { 5926322f4bdSBarry Smith for (i=istart; i<i2end; i++) { 5936322f4bdSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 594b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 595b2f1ab58SBarry Smith i1++; 5966322f4bdSBarry Smith } else if (i2<i2end) { 597b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 598b2f1ab58SBarry Smith i2++; 5996322f4bdSBarry Smith } else { 600b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 601b2f1ab58SBarry Smith i1++; 602b2f1ab58SBarry Smith } 603b2f1ab58SBarry Smith } 604b2f1ab58SBarry Smith } 605b2f1ab58SBarry Smith } 606b2f1ab58SBarry Smith 607c328eeadSBarry Smith /* Swap the two array sets */ 608b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 609b2f1ab58SBarry Smith vswap = vhlp2; vhlp2 = vhlp1; vhlp1 = vswap; 610b2f1ab58SBarry Smith } 611b2f1ab58SBarry Smith 612c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 6136322f4bdSBarry Smith if (ihlp2 != icol) { 6142205254eSKarl Rupp for (i=0; i<nnz; i++) icol[i] = ihlp2[i]; 6156322f4bdSBarry Smith if (val) { 6162205254eSKarl Rupp for (i=0; i<nnz; i++) val[i] = vhlp2[i]; 617b2f1ab58SBarry Smith } 618b2f1ab58SBarry Smith } 619b2f1ab58SBarry Smith 6209566063dSJacob Faibussowitsch PetscCall(PetscFree(ialloc)); 6219566063dSJacob Faibussowitsch if (val) PetscCall(PetscFree(valloc)); 622b2f1ab58SBarry Smith PetscFunctionReturn(0); 623b2f1ab58SBarry Smith } 624b2f1ab58SBarry Smith 625b2f1ab58SBarry Smith /* 626b2f1ab58SBarry Smith spbas_apply_reordering_rows: 627b2f1ab58SBarry Smith apply the given reordering to the rows: matrix_A = matrix_A(perm,:); 628b2f1ab58SBarry Smith */ 629b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_rows(spbas_matrix *matrix_A, const PetscInt *permutation) 630b2f1ab58SBarry Smith { 631b2f1ab58SBarry Smith PetscInt i,j,ip; 632b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 633b2f1ab58SBarry Smith PetscInt * row_nnz; 634b2f1ab58SBarry Smith PetscInt **icols; 635ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 6360298fd71SBarry Smith PetscScalar **vals = NULL; 637b2f1ab58SBarry Smith 638b2f1ab58SBarry Smith PetscFunctionBegin; 639*08401ef6SPierre Jolivet PetscCheck(matrix_A->col_idx_type == SPBAS_DIAGONAL_OFFSETS,PETSC_COMM_SELF, PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern"); 640b2f1ab58SBarry Smith 6416322f4bdSBarry Smith if (do_values) { 6429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows, &vals)); 643b2f1ab58SBarry Smith } 6449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows, &row_nnz)); 6459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows, &icols)); 646b2f1ab58SBarry Smith 6476322f4bdSBarry Smith for (i=0; i<nrows; i++) { 648b2f1ab58SBarry Smith ip = permutation[i]; 6492205254eSKarl Rupp if (do_values) vals[i] = matrix_A->values[ip]; 650b2f1ab58SBarry Smith icols[i] = matrix_A->icols[ip]; 651b2f1ab58SBarry Smith row_nnz[i] = matrix_A->row_nnz[ip]; 6522205254eSKarl Rupp for (j=0; j<row_nnz[i]; j++) icols[i][j] += ip-i; 653b2f1ab58SBarry Smith } 654b2f1ab58SBarry Smith 6559566063dSJacob Faibussowitsch if (do_values) PetscCall(PetscFree(matrix_A->values)); 6569566063dSJacob Faibussowitsch PetscCall(PetscFree(matrix_A->icols)); 6579566063dSJacob Faibussowitsch PetscCall(PetscFree(matrix_A->row_nnz)); 658b2f1ab58SBarry Smith 6592205254eSKarl Rupp if (do_values) matrix_A->values = vals; 660b2f1ab58SBarry Smith matrix_A->icols = icols; 661b2f1ab58SBarry Smith matrix_A->row_nnz = row_nnz; 662b2f1ab58SBarry Smith PetscFunctionReturn(0); 663b2f1ab58SBarry Smith } 664b2f1ab58SBarry Smith 665b2f1ab58SBarry Smith /* 666b2f1ab58SBarry Smith spbas_apply_reordering_cols: 667b2f1ab58SBarry Smith apply the given reordering to the columns: matrix_A(:,perm) = matrix_A; 668b2f1ab58SBarry Smith */ 669b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_cols(spbas_matrix *matrix_A,const PetscInt *permutation) 670b2f1ab58SBarry Smith { 671b2f1ab58SBarry Smith PetscInt i,j; 672b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 673b2f1ab58SBarry Smith PetscInt row_nnz; 674b2f1ab58SBarry Smith PetscInt *icols; 675ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 6760298fd71SBarry Smith PetscScalar *vals = NULL; 677b2f1ab58SBarry Smith 678b2f1ab58SBarry Smith PetscFunctionBegin; 679*08401ef6SPierre Jolivet PetscCheck(matrix_A->col_idx_type == SPBAS_DIAGONAL_OFFSETS,PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "must have diagonal offsets in pattern"); 680b2f1ab58SBarry Smith 6816322f4bdSBarry Smith for (i=0; i<nrows; i++) { 682b2f1ab58SBarry Smith icols = matrix_A->icols[i]; 683b2f1ab58SBarry Smith row_nnz = matrix_A->row_nnz[i]; 6842205254eSKarl Rupp if (do_values) vals = matrix_A->values[i]; 685b2f1ab58SBarry Smith 6866322f4bdSBarry Smith for (j=0; j<row_nnz; j++) { 687b2f1ab58SBarry Smith icols[j] = permutation[i+icols[j]]-i; 688b2f1ab58SBarry Smith } 6899566063dSJacob Faibussowitsch PetscCall(spbas_mergesort(row_nnz, icols, vals)); 690b2f1ab58SBarry Smith } 691b2f1ab58SBarry Smith PetscFunctionReturn(0); 692b2f1ab58SBarry Smith } 693b2f1ab58SBarry Smith 694b2f1ab58SBarry Smith /* 695b2f1ab58SBarry Smith spbas_apply_reordering: 696b2f1ab58SBarry Smith apply the given reordering: matrix_A(perm,perm) = matrix_A; 697b2f1ab58SBarry Smith */ 698b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering(spbas_matrix *matrix_A, const PetscInt *permutation, const PetscInt * inv_perm) 699b2f1ab58SBarry Smith { 700b2f1ab58SBarry Smith PetscFunctionBegin; 7019566063dSJacob Faibussowitsch PetscCall(spbas_apply_reordering_rows(matrix_A, inv_perm)); 7029566063dSJacob Faibussowitsch PetscCall(spbas_apply_reordering_cols(matrix_A, permutation)); 703b2f1ab58SBarry Smith PetscFunctionReturn(0); 704b2f1ab58SBarry Smith } 705b2f1ab58SBarry Smith 706b2f1ab58SBarry Smith PetscErrorCode spbas_pattern_only(PetscInt nrows, PetscInt ncols, PetscInt *ai, PetscInt *aj, spbas_matrix * result) 707b2f1ab58SBarry Smith { 708b2f1ab58SBarry Smith spbas_matrix retval; 709b2f1ab58SBarry Smith PetscInt i, j, i0, r_nnz; 710b2f1ab58SBarry Smith 711b2f1ab58SBarry Smith PetscFunctionBegin; 712c328eeadSBarry Smith /* Copy input values */ 713b2f1ab58SBarry Smith retval.nrows = nrows; 714b2f1ab58SBarry Smith retval.ncols = ncols; 715b2f1ab58SBarry Smith retval.nnz = ai[nrows]; 716b2f1ab58SBarry Smith 717b2f1ab58SBarry Smith retval.block_data = PETSC_TRUE; 718b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 719b2f1ab58SBarry Smith 720c328eeadSBarry Smith /* Allocate output matrix */ 7219566063dSJacob Faibussowitsch PetscCall(spbas_allocate_pattern(&retval, PETSC_FALSE)); 7222205254eSKarl Rupp for (i=0; i<nrows; i++) retval.row_nnz[i] = ai[i+1]-ai[i]; 7239566063dSJacob Faibussowitsch PetscCall(spbas_allocate_data(&retval)); 724c328eeadSBarry Smith /* Copy the structure */ 7256322f4bdSBarry Smith for (i = 0; i<retval.nrows; i++) { 726b2f1ab58SBarry Smith i0 = ai[i]; 727b2f1ab58SBarry Smith r_nnz = ai[i+1]-i0; 728b2f1ab58SBarry Smith 7296322f4bdSBarry Smith for (j=0; j<r_nnz; j++) { 730b2f1ab58SBarry Smith retval.icols[i][j] = aj[i0+j]-i; 731b2f1ab58SBarry Smith } 732b2f1ab58SBarry Smith } 733b2f1ab58SBarry Smith *result = retval; 734b2f1ab58SBarry Smith PetscFunctionReturn(0); 735b2f1ab58SBarry Smith } 736b2f1ab58SBarry Smith 737b2f1ab58SBarry Smith /* 738b2f1ab58SBarry Smith spbas_mark_row_power: 739b2f1ab58SBarry Smith Mark the columns in row 'row' which are nonzero in 740b2f1ab58SBarry Smith matrix^2log(marker). 741b2f1ab58SBarry Smith */ 742be332245SKarl Rupp PetscErrorCode spbas_mark_row_power(PetscInt *iwork, /* marker-vector */ 743c328eeadSBarry Smith PetscInt row, /* row for which the columns are marked */ 744c328eeadSBarry Smith spbas_matrix * in_matrix, /* matrix for which the power is being calculated */ 745c328eeadSBarry Smith PetscInt marker, /* marker-value: 2^power */ 746c328eeadSBarry Smith PetscInt minmrk, /* lower bound for marked points */ 747c328eeadSBarry Smith PetscInt maxmrk) /* upper bound for marked points */ 748b2f1ab58SBarry Smith { 749b2f1ab58SBarry Smith PetscInt i,j, nnz; 750b2f1ab58SBarry Smith 751b2f1ab58SBarry Smith PetscFunctionBegin; 752b2f1ab58SBarry Smith nnz = in_matrix->row_nnz[row]; 753b2f1ab58SBarry Smith 754c328eeadSBarry Smith /* For higher powers, call this function recursively */ 7556322f4bdSBarry Smith if (marker>1) { 7566322f4bdSBarry Smith for (i=0; i<nnz; i++) { 757b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 7586322f4bdSBarry Smith if (minmrk<=j && j<maxmrk && iwork[j] < marker) { 7599566063dSJacob Faibussowitsch PetscCall(spbas_mark_row_power(iwork, row + in_matrix->icols[row][i],in_matrix, marker/2,minmrk,maxmrk)); 760b2f1ab58SBarry Smith iwork[j] |= marker; 761b2f1ab58SBarry Smith } 762b2f1ab58SBarry Smith } 7636322f4bdSBarry Smith } else { 764c328eeadSBarry Smith /* Mark the columns reached */ 7656322f4bdSBarry Smith for (i=0; i<nnz; i++) { 766b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 7672205254eSKarl Rupp if (minmrk<=j && j<maxmrk) iwork[j] |= 1; 768b2f1ab58SBarry Smith } 769b2f1ab58SBarry Smith } 770b2f1ab58SBarry Smith PetscFunctionReturn(0); 771b2f1ab58SBarry Smith } 772b2f1ab58SBarry Smith 773b2f1ab58SBarry Smith /* 774b2f1ab58SBarry Smith spbas_power 775b2f1ab58SBarry Smith Calculate sparseness patterns for incomplete Cholesky decompositions 776b2f1ab58SBarry Smith of a given order: (almost) all nonzeros of the matrix^(order+1) which 777b2f1ab58SBarry Smith are inside the band width are found and stored in the output sparseness 778b2f1ab58SBarry Smith pattern. 779b2f1ab58SBarry Smith */ 780b2f1ab58SBarry Smith PetscErrorCode spbas_power(spbas_matrix in_matrix,PetscInt power, spbas_matrix * result) 781b2f1ab58SBarry Smith { 782b2f1ab58SBarry Smith spbas_matrix retval; 783b2f1ab58SBarry Smith PetscInt nrows = in_matrix.nrows; 784b2f1ab58SBarry Smith PetscInt ncols = in_matrix.ncols; 785b2f1ab58SBarry Smith PetscInt i, j, kend; 786b2f1ab58SBarry Smith PetscInt nnz, inz; 787b2f1ab58SBarry Smith PetscInt *iwork; 788b2f1ab58SBarry Smith PetscInt marker; 789b2f1ab58SBarry Smith PetscInt maxmrk=0; 790b2f1ab58SBarry Smith 791b2f1ab58SBarry Smith PetscFunctionBegin; 792*08401ef6SPierre Jolivet PetscCheck(in_matrix.col_idx_type == SPBAS_DIAGONAL_OFFSETS,PETSC_COMM_SELF, PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern"); 793*08401ef6SPierre Jolivet PetscCheck(ncols == nrows,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Dimension error"); 7942c71b3e2SJacob Faibussowitsch PetscCheckFalse(in_matrix.values,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input array must be sparseness pattern (no values)"); 795*08401ef6SPierre Jolivet PetscCheck(power>0,PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Power must be 1 or up"); 796b2f1ab58SBarry Smith 797c328eeadSBarry Smith /* Copy input values*/ 798b2f1ab58SBarry Smith retval.nrows = ncols; 799b2f1ab58SBarry Smith retval.ncols = nrows; 800b2f1ab58SBarry Smith retval.nnz = 0; 801b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 802b2f1ab58SBarry Smith retval.block_data = PETSC_FALSE; 803b2f1ab58SBarry Smith 804c328eeadSBarry Smith /* Allocate sparseness pattern */ 8059566063dSJacob Faibussowitsch PetscCall(spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE)); 806b2f1ab58SBarry Smith 807580bdb30SBarry Smith /* Allocate marker array: note sure the max needed so use the max of the two */ 8089566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(PetscMax(ncols,nrows), &iwork)); 809b2f1ab58SBarry Smith 810c328eeadSBarry Smith /* Calculate marker values */ 8112205254eSKarl Rupp marker = 1; for (i=1; i<power; i++) marker*=2; 812b2f1ab58SBarry Smith 8136322f4bdSBarry Smith for (i=0; i<nrows; i++) { 814c328eeadSBarry Smith /* Calculate the pattern for each row */ 815b2f1ab58SBarry Smith 816b2f1ab58SBarry Smith nnz = in_matrix.row_nnz[i]; 817b2f1ab58SBarry Smith kend = i+in_matrix.icols[i][nnz-1]; 8182205254eSKarl Rupp if (maxmrk<=kend) maxmrk=kend+1; 8199566063dSJacob Faibussowitsch PetscCall(spbas_mark_row_power(iwork, i, &in_matrix, marker, i, maxmrk)); 820b2f1ab58SBarry Smith 821c328eeadSBarry Smith /* Count the columns*/ 822b2f1ab58SBarry Smith nnz = 0; 8232205254eSKarl Rupp for (j=i; j<maxmrk; j++) nnz+= (iwork[j]!=0); 824b2f1ab58SBarry Smith 825c328eeadSBarry Smith /* Allocate the column indices */ 826b2f1ab58SBarry Smith retval.row_nnz[i] = nnz; 8279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz,&retval.icols[i])); 828b2f1ab58SBarry Smith 829c328eeadSBarry Smith /* Administrate the column indices */ 830b2f1ab58SBarry Smith inz = 0; 8316322f4bdSBarry Smith for (j=i; j<maxmrk; j++) { 8326322f4bdSBarry Smith if (iwork[j]) { 833b2f1ab58SBarry Smith retval.icols[i][inz] = j-i; 834b2f1ab58SBarry Smith inz++; 835b2f1ab58SBarry Smith iwork[j] = 0; 836b2f1ab58SBarry Smith } 837b2f1ab58SBarry Smith } 838b2f1ab58SBarry Smith retval.nnz += nnz; 839b2f1ab58SBarry Smith }; 8409566063dSJacob Faibussowitsch PetscCall(PetscFree(iwork)); 841b2f1ab58SBarry Smith *result = retval; 842b2f1ab58SBarry Smith PetscFunctionReturn(0); 843b2f1ab58SBarry Smith } 844b2f1ab58SBarry Smith 845b2f1ab58SBarry Smith /* 846b2f1ab58SBarry Smith spbas_keep_upper: 847b2f1ab58SBarry Smith remove the lower part of the matrix: keep the upper part 848b2f1ab58SBarry Smith */ 849b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper(spbas_matrix * inout_matrix) 850b2f1ab58SBarry Smith { 851b2f1ab58SBarry Smith PetscInt i, j; 852b2f1ab58SBarry Smith PetscInt jstart; 853b2f1ab58SBarry Smith 854b2f1ab58SBarry Smith PetscFunctionBegin; 8555f80ce2aSJacob Faibussowitsch PetscCheck(!inout_matrix->block_data,PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Not yet for block data matrices"); 8566322f4bdSBarry Smith for (i=0; i<inout_matrix->nrows; i++) { 8576322f4bdSBarry Smith for (jstart=0; (jstart<inout_matrix->row_nnz[i]) && (inout_matrix->icols[i][jstart]<0); jstart++) {} 8586322f4bdSBarry Smith if (jstart>0) { 8596322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 8606322f4bdSBarry Smith inout_matrix->icols[i][j] = inout_matrix->icols[i][j+jstart]; 861b2f1ab58SBarry Smith } 862b2f1ab58SBarry Smith 8636322f4bdSBarry Smith if (inout_matrix->values) { 8646322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 8656322f4bdSBarry Smith inout_matrix->values[i][j] = inout_matrix->values[i][j+jstart]; 866b2f1ab58SBarry Smith } 867b2f1ab58SBarry Smith } 868b2f1ab58SBarry Smith 869b2f1ab58SBarry Smith inout_matrix->row_nnz[i] -= jstart; 870b2f1ab58SBarry Smith 8716322f4bdSBarry Smith inout_matrix->icols[i] = (PetscInt*) realloc((void*) inout_matrix->icols[i], inout_matrix->row_nnz[i]*sizeof(PetscInt)); 872b2f1ab58SBarry Smith 8736322f4bdSBarry Smith if (inout_matrix->values) { 8746322f4bdSBarry Smith inout_matrix->values[i] = (PetscScalar*) realloc((void*) inout_matrix->values[i], inout_matrix->row_nnz[i]*sizeof(PetscScalar)); 875b2f1ab58SBarry Smith } 876b2f1ab58SBarry Smith inout_matrix->nnz -= jstart; 877b2f1ab58SBarry Smith } 878b2f1ab58SBarry Smith } 879b2f1ab58SBarry Smith PetscFunctionReturn(0); 880b2f1ab58SBarry Smith } 881