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: 10845006b9SBarry Smith + -pc_factor_levels <l> 11b7c853c4SBarry Smith - -pc_factor_drop_tolerance 12845006b9SBarry Smith 13845006b9SBarry Smith Level: beginner 14845006b9SBarry Smith 15845006b9SBarry Smith Contributed by: Bas van 't Hof 16845006b9SBarry Smith 17b7c853c4SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, PCFactorSetLevels(), PCFactorSetDropTolerance() 18845006b9SBarry Smith 19845006b9SBarry Smith M*/ 209ccc4a9bSBarry Smith 21b2f1ab58SBarry Smith /* 22b2f1ab58SBarry Smith spbas_memory_requirement: 23b2f1ab58SBarry Smith Calculate the number of bytes needed to store tha matrix 24b2f1ab58SBarry Smith */ 25b2f1ab58SBarry Smith #undef __FUNCT__ 26b2f1ab58SBarry Smith #define __FUNCT__ "spbas_memory_requirement" 27b2f1ab58SBarry Smith long int spbas_memory_requirement(spbas_matrix matrix) 28b2f1ab58SBarry Smith { 29*2205254eSKarl Rupp long int memreq = 6 * sizeof(PetscInt) + /* nrows, ncols, nnz, n_alloc_icol, n_alloc_val, col_idx_type */ 30ace3abfcSBarry Smith sizeof(PetscBool) + /* block_data */ 31c328eeadSBarry Smith sizeof(PetscScalar**) + /* values */ 32c328eeadSBarry Smith sizeof(PetscScalar*) + /* alloc_val */ 33c328eeadSBarry Smith 2 * sizeof(int**) + /* icols, icols0 */ 34c328eeadSBarry Smith 2 * sizeof(PetscInt*) + /* row_nnz, alloc_icol */ 35c328eeadSBarry Smith matrix.nrows * sizeof(PetscInt) + /* row_nnz[*] */ 36c328eeadSBarry Smith matrix.nrows * sizeof(PetscInt*); /* icols[*] */ 37b2f1ab58SBarry Smith 38c328eeadSBarry Smith /* icol0[*] */ 39*2205254eSKarl Rupp if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) memreq += matrix.nrows * sizeof(PetscInt); 40b2f1ab58SBarry Smith 41c328eeadSBarry Smith /* icols[*][*] */ 42*2205254eSKarl Rupp if (matrix.block_data) memreq += matrix.n_alloc_icol * sizeof(PetscInt); 43*2205254eSKarl Rupp else memreq += matrix.nnz * sizeof(PetscInt); 44b2f1ab58SBarry Smith 454e1ff37aSBarry Smith if (matrix.values) { 46c328eeadSBarry Smith memreq += matrix.nrows * sizeof(PetscScalar*); /* values[*] */ 47c328eeadSBarry Smith /* values[*][*] */ 48*2205254eSKarl Rupp if (matrix.block_data) memreq += matrix.n_alloc_val * sizeof(PetscScalar); 49*2205254eSKarl Rupp else memreq += matrix.nnz * sizeof(PetscScalar); 50b2f1ab58SBarry Smith } 51b2f1ab58SBarry Smith return memreq; 52b2f1ab58SBarry Smith } 53b2f1ab58SBarry Smith 54b2f1ab58SBarry Smith /* 55b2f1ab58SBarry Smith spbas_allocate_pattern: 56b2f1ab58SBarry Smith allocate the pattern arrays row_nnz, icols and optionally values 57b2f1ab58SBarry Smith */ 58b2f1ab58SBarry Smith #undef __FUNCT__ 59b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_pattern" 60ace3abfcSBarry Smith PetscErrorCode spbas_allocate_pattern(spbas_matrix * result, PetscBool do_values) 61b2f1ab58SBarry Smith { 62b2f1ab58SBarry Smith PetscErrorCode ierr; 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 */ 68b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->row_nnz);CHKERRQ(ierr); 69b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt*),&result->icols);CHKERRQ(ierr); 70b2f1ab58SBarry Smith 71c328eeadSBarry Smith /* If offsets are given wrt an array, create array */ 724e1ff37aSBarry Smith if (col_idx_type == SPBAS_OFFSET_ARRAY) { 73b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->icol0);CHKERRQ(ierr); 744e1ff37aSBarry Smith } else { 75c328eeadSBarry Smith result->icol0 = PETSC_NULL; 76b2f1ab58SBarry Smith } 77b2f1ab58SBarry Smith 78c328eeadSBarry Smith /* If values are given, allocate values array */ 794e1ff37aSBarry Smith if (do_values) { 80c328eeadSBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscScalar*),&result->values);CHKERRQ(ierr); 814e1ff37aSBarry Smith } else { 82c328eeadSBarry Smith result->values = PETSC_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 #undef __FUNCT__ 96b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_data" 97b2f1ab58SBarry Smith PetscErrorCode spbas_allocate_data(spbas_matrix * result) 98b2f1ab58SBarry Smith { 99b2f1ab58SBarry Smith PetscInt i; 100b2f1ab58SBarry Smith PetscInt nnz = result->nnz; 101b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 102b2f1ab58SBarry Smith PetscInt r_nnz; 103b2f1ab58SBarry Smith PetscErrorCode ierr; 104ace3abfcSBarry Smith PetscBool do_values = (result->values != PETSC_NULL) ? PETSC_TRUE : PETSC_FALSE; 105ace3abfcSBarry Smith PetscBool block_data = result->block_data; 106b2f1ab58SBarry Smith 107b2f1ab58SBarry Smith PetscFunctionBegin; 1084e1ff37aSBarry Smith if (block_data) { 109c328eeadSBarry Smith /* Allocate the column number array and point to it */ 110b2f1ab58SBarry Smith result->n_alloc_icol = nnz; 111*2205254eSKarl Rupp 112c328eeadSBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt), &result->alloc_icol);CHKERRQ(ierr); 113*2205254eSKarl Rupp 114b2f1ab58SBarry Smith result->icols[0] = result->alloc_icol; 1154e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 116b2f1ab58SBarry Smith result->icols[i] = result->icols[i-1] + result->row_nnz[i-1]; 117b2f1ab58SBarry Smith } 118b2f1ab58SBarry Smith 119c328eeadSBarry Smith /* Allocate the value array and point to it */ 1204e1ff37aSBarry Smith if (do_values) { 121b2f1ab58SBarry Smith result->n_alloc_val = nnz; 122*2205254eSKarl Rupp 123c328eeadSBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscScalar), &result->alloc_val);CHKERRQ(ierr); 124*2205254eSKarl Rupp 125b2f1ab58SBarry Smith result->values[0] = result->alloc_val; 1264e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 127b2f1ab58SBarry Smith result->values[i] = result->values[i-1] + result->row_nnz[i-1]; 128b2f1ab58SBarry Smith } 129b2f1ab58SBarry Smith } 1304e1ff37aSBarry Smith } else { 1314e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 132b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 133c328eeadSBarry Smith ierr = PetscMalloc(r_nnz*sizeof(PetscInt), &result->icols[i]);CHKERRQ(ierr); 134b2f1ab58SBarry Smith } 1354e1ff37aSBarry Smith if (do_values) { 1364e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 137b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 1384e1ff37aSBarry Smith ierr = PetscMalloc(r_nnz*sizeof(PetscScalar), &result->values[i]);CHKERRQ(ierr); 139b2f1ab58SBarry Smith } 140b2f1ab58SBarry Smith } 141b2f1ab58SBarry Smith } 142b2f1ab58SBarry Smith PetscFunctionReturn(0); 143b2f1ab58SBarry Smith } 144b2f1ab58SBarry Smith 145b2f1ab58SBarry Smith /* 146b2f1ab58SBarry Smith spbas_row_order_icol 147b2f1ab58SBarry Smith determine if row i1 should come 148b2f1ab58SBarry Smith + before row i2 in the sorted rows (return -1), 149b2f1ab58SBarry Smith + after (return 1) 150b2f1ab58SBarry Smith + is identical (return 0). 151b2f1ab58SBarry Smith */ 152b2f1ab58SBarry Smith #undef __FUNCT__ 153b2f1ab58SBarry Smith #define __FUNCT__ "spbas_row_order_icol" 154b2f1ab58SBarry Smith int spbas_row_order_icol(PetscInt i1, PetscInt i2, PetscInt *irow_in, PetscInt *icol_in,PetscInt col_idx_type) 155b2f1ab58SBarry Smith { 156b2f1ab58SBarry Smith PetscInt j; 157b2f1ab58SBarry Smith PetscInt nnz1 = irow_in[i1+1] - irow_in[i1]; 158b2f1ab58SBarry Smith PetscInt nnz2 = irow_in[i2+1] - irow_in[i2]; 159b2f1ab58SBarry Smith PetscInt * icol1 = &icol_in[irow_in[i1]]; 160b2f1ab58SBarry Smith PetscInt * icol2 = &icol_in[irow_in[i2]]; 161b2f1ab58SBarry Smith 162*2205254eSKarl Rupp if (nnz1<nnz2) return -1; 163*2205254eSKarl Rupp if (nnz1>nnz2) return 1; 164b2f1ab58SBarry Smith 1654e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 1664e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 167*2205254eSKarl Rupp if (icol1[j]< icol2[j]) return -1; 168*2205254eSKarl Rupp if (icol1[j]> icol2[j]) return 1; 169b2f1ab58SBarry Smith } 1704e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 1714e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 172*2205254eSKarl Rupp if (icol1[j]-i1< icol2[j]-i2) return -1; 173*2205254eSKarl Rupp if (icol1[j]-i1> icol2[j]-i2) return 1; 174b2f1ab58SBarry Smith } 1754e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 1764e1ff37aSBarry Smith for (j=1; j<nnz1; j++) { 177*2205254eSKarl Rupp if (icol1[j]-icol1[0] < icol2[j]-icol2[0]) return -1; 178*2205254eSKarl Rupp if (icol1[j]-icol1[0] > icol2[j]-icol2[0]) return 1; 179b2f1ab58SBarry Smith } 180b2f1ab58SBarry Smith } 181b2f1ab58SBarry Smith return 0; 182b2f1ab58SBarry Smith } 183b2f1ab58SBarry Smith 184b2f1ab58SBarry Smith /* 185b2f1ab58SBarry Smith spbas_mergesort_icols: 186b2f1ab58SBarry Smith return a sorting of the rows in which identical sparseness patterns are 187b2f1ab58SBarry Smith next to each other 188b2f1ab58SBarry Smith */ 189b2f1ab58SBarry Smith #undef __FUNCT__ 190b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort_icols" 191b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort_icols(PetscInt nrows, PetscInt * irow_in, PetscInt * icol_in,PetscInt col_idx_type, PetscInt *isort) 192b2f1ab58SBarry Smith { 193b2f1ab58SBarry Smith PetscErrorCode ierr; 194c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 195c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 196c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 197c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 198c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 199c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 200c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 201b2f1ab58SBarry Smith 202b2f1ab58SBarry Smith PetscFunctionBegin; 203b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&ialloc);CHKERRQ(ierr); 204b2f1ab58SBarry Smith 205b2f1ab58SBarry Smith ihlp1 = ialloc; 206b2f1ab58SBarry Smith ihlp2 = isort; 207b2f1ab58SBarry Smith 208c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 2094e1ff37aSBarry Smith for (istep=1; istep<nrows; istep*=2) { 210c328eeadSBarry Smith /* 211c328eeadSBarry Smith Combine sorted parts 212c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 213c328eeadSBarry Smith of ihlp2 and vhlp2 214c328eeadSBarry Smith 215c328eeadSBarry Smith into one sorted part 216c328eeadSBarry Smith istart:istart+2*istep-1 217c328eeadSBarry Smith of ihlp1 and vhlp1 218c328eeadSBarry Smith */ 2194e1ff37aSBarry Smith for (istart=0; istart<nrows; istart+=2*istep) { 220c328eeadSBarry Smith /* Set counters and bound array part endings */ 221*2205254eSKarl Rupp i1=istart; i1end = i1+istep; if (i1end>nrows) i1end=nrows; 222*2205254eSKarl Rupp i2=istart+istep; i2end = i2+istep; if (i2end>nrows) i2end=nrows; 223b2f1ab58SBarry Smith 224c328eeadSBarry Smith /* Merge the two array parts */ 2254e1ff37aSBarry Smith for (i=istart; i<i2end; i++) { 2264e1ff37aSBarry Smith if (i1<i1end && i2<i2end && spbas_row_order_icol(ihlp2[i1], ihlp2[i2], irow_in, icol_in, col_idx_type) < 0) { 227b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 228b2f1ab58SBarry Smith i1++; 2294e1ff37aSBarry Smith } else if (i2<i2end) { 230b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 231b2f1ab58SBarry Smith i2++; 2324e1ff37aSBarry Smith } else { 233b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 234b2f1ab58SBarry Smith i1++; 235b2f1ab58SBarry Smith } 236b2f1ab58SBarry Smith } 237b2f1ab58SBarry Smith } 238b2f1ab58SBarry Smith 239c328eeadSBarry Smith /* Swap the two array sets */ 240b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 241b2f1ab58SBarry Smith } 242b2f1ab58SBarry Smith 243c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 2444e1ff37aSBarry Smith if (ihlp2 != isort) { 245*2205254eSKarl Rupp for (i=0; i<nrows; i++) isort[i] = ihlp2[i]; 246b2f1ab58SBarry Smith } 247b2f1ab58SBarry Smith ierr = PetscFree(ialloc);CHKERRQ(ierr); 248b2f1ab58SBarry Smith PetscFunctionReturn(0); 249b2f1ab58SBarry Smith } 250b2f1ab58SBarry Smith 251b2f1ab58SBarry Smith 252b2f1ab58SBarry Smith 253b2f1ab58SBarry Smith /* 254b2f1ab58SBarry Smith spbas_compress_pattern: 255b2f1ab58SBarry Smith calculate a compressed sparseness pattern for a sparseness pattern 256b2f1ab58SBarry Smith given in compressed row storage. The compressed sparseness pattern may 257b2f1ab58SBarry Smith require (much) less memory. 258b2f1ab58SBarry Smith */ 259b2f1ab58SBarry Smith #undef __FUNCT__ 260b2f1ab58SBarry Smith #define __FUNCT__ "spbas_compress_pattern" 2614e1ff37aSBarry 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) 262b2f1ab58SBarry Smith { 263b2f1ab58SBarry Smith PetscInt nnz = irow_in[nrows]; 264b2f1ab58SBarry Smith long int mem_orig = (nrows + nnz) * sizeof(PetscInt); 265b2f1ab58SBarry Smith long int mem_compressed; 266b2f1ab58SBarry Smith PetscErrorCode ierr; 267b2f1ab58SBarry Smith PetscInt *isort; 268b2f1ab58SBarry Smith PetscInt *icols; 269b2f1ab58SBarry Smith PetscInt row_nnz; 270b2f1ab58SBarry Smith PetscInt *ipoint; 271ace3abfcSBarry Smith PetscBool *used; 272b2f1ab58SBarry Smith PetscInt ptr; 273b2f1ab58SBarry Smith PetscInt i,j; 274ace3abfcSBarry Smith const PetscBool no_values = PETSC_FALSE; 275b2f1ab58SBarry Smith 276b2f1ab58SBarry Smith PetscFunctionBegin; 277c328eeadSBarry Smith /* Allocate the structure of the new matrix */ 278b2f1ab58SBarry Smith B->nrows = nrows; 279b2f1ab58SBarry Smith B->ncols = ncols; 280b2f1ab58SBarry Smith B->nnz = nnz; 281b2f1ab58SBarry Smith B->col_idx_type = col_idx_type; 282b2f1ab58SBarry Smith B->block_data = PETSC_TRUE; 283*2205254eSKarl Rupp 284b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(B, no_values);CHKERRQ(ierr); 285b2f1ab58SBarry Smith 286c328eeadSBarry Smith /* When using an offset array, set it */ 2874e1ff37aSBarry Smith if (col_idx_type==SPBAS_OFFSET_ARRAY) { 288*2205254eSKarl Rupp for (i=0; i<nrows; i++) B->icol0[i] = icol_in[irow_in[i]]; 289b2f1ab58SBarry Smith } 290b2f1ab58SBarry Smith 291c328eeadSBarry Smith /* Allocate the ordering for the rows */ 292b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&isort);CHKERRQ(ierr); 293b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&ipoint);CHKERRQ(ierr); 294ace3abfcSBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscBool),&used);CHKERRQ(ierr); 295b2f1ab58SBarry Smith 296c328eeadSBarry Smith /* Initialize the sorting */ 297ace3abfcSBarry Smith ierr = PetscMemzero((void*) used, nrows*sizeof(PetscBool));CHKERRQ(ierr); 2984e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 299b2f1ab58SBarry Smith B->row_nnz[i] = irow_in[i+1]-irow_in[i]; 300b2f1ab58SBarry Smith isort[i] = i; 301b2f1ab58SBarry Smith ipoint[i] = i; 302b2f1ab58SBarry Smith } 303b2f1ab58SBarry Smith 304c328eeadSBarry Smith /* Sort the rows so that identical columns will be next to each other */ 305b2f1ab58SBarry Smith ierr = spbas_mergesort_icols(nrows, irow_in, icol_in, col_idx_type, isort);CHKERRQ(ierr); 30671d55d18SBarry Smith ierr = PetscInfo(PETSC_NULL,"Rows have been sorted for patterns\n");CHKERRQ(ierr); 307b2f1ab58SBarry Smith 308c328eeadSBarry Smith /* Replace identical rows with the first one in the list */ 3094e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 3104e1ff37aSBarry Smith if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) { 311b2f1ab58SBarry Smith ipoint[isort[i]] = ipoint[isort[i-1]]; 312b2f1ab58SBarry Smith } 313b2f1ab58SBarry Smith } 314b2f1ab58SBarry Smith 315c328eeadSBarry Smith /* Collect the rows which are used*/ 316*2205254eSKarl Rupp for (i=0; i<nrows; i++) used[ipoint[i]] = PETSC_TRUE; 317b2f1ab58SBarry Smith 318c328eeadSBarry Smith /* Calculate needed memory */ 319b2f1ab58SBarry Smith B->n_alloc_icol = 0; 3204e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 321*2205254eSKarl Rupp if (used[i]) B->n_alloc_icol += B->row_nnz[i]; 322b2f1ab58SBarry Smith } 32371d55d18SBarry Smith ierr = PetscMalloc(B->n_alloc_icol*sizeof(PetscInt),&B->alloc_icol);CHKERRQ(ierr); 324b2f1ab58SBarry Smith 325c328eeadSBarry Smith /* Fill in the diagonal offsets for the rows which store their own data */ 326b2f1ab58SBarry Smith ptr = 0; 3274e1ff37aSBarry Smith for (i=0; i<B->nrows; i++) { 3284e1ff37aSBarry Smith if (used[i]) { 329b2f1ab58SBarry Smith B->icols[i] = &B->alloc_icol[ptr]; 330b2f1ab58SBarry Smith icols = &icol_in[irow_in[i]]; 331b2f1ab58SBarry Smith row_nnz = B->row_nnz[i]; 3324e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 3334e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 334b2f1ab58SBarry Smith B->icols[i][j] = icols[j]; 335b2f1ab58SBarry Smith } 3364e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 3374e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 338b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-i; 339b2f1ab58SBarry Smith } 3404e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 3414e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 342b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-icols[0]; 343b2f1ab58SBarry Smith } 344b2f1ab58SBarry Smith } 345b2f1ab58SBarry Smith ptr += B->row_nnz[i]; 346b2f1ab58SBarry Smith } 347b2f1ab58SBarry Smith } 348b2f1ab58SBarry Smith 349c328eeadSBarry Smith /* Point to the right places for all data */ 3504e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 351b2f1ab58SBarry Smith B->icols[i] = B->icols[ipoint[i]]; 352b2f1ab58SBarry Smith } 353c328eeadSBarry Smith ierr = PetscInfo(PETSC_NULL,"Row patterns have been compressed\n");CHKERRQ(ierr); 354c328eeadSBarry Smith ierr = PetscInfo1(PETSC_NULL," (%G nonzeros per row)\n", (PetscReal) nnz / (PetscReal) nrows);CHKERRQ(ierr); 355b2f1ab58SBarry Smith 356b2f1ab58SBarry Smith ierr=PetscFree(isort);CHKERRQ(ierr); 357b2f1ab58SBarry Smith ierr=PetscFree(used);CHKERRQ(ierr); 358b2f1ab58SBarry Smith ierr=PetscFree(ipoint);CHKERRQ(ierr); 359b2f1ab58SBarry Smith 360b2f1ab58SBarry Smith mem_compressed = spbas_memory_requirement(*B); 3614e1ff37aSBarry Smith *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig; 362b2f1ab58SBarry Smith PetscFunctionReturn(0); 363b2f1ab58SBarry Smith } 364b2f1ab58SBarry Smith 365b2f1ab58SBarry Smith /* 366b2f1ab58SBarry Smith spbas_incomplete_cholesky 367b2f1ab58SBarry Smith Incomplete Cholesky decomposition 368b2f1ab58SBarry Smith */ 369c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas_cholesky.h> 370b2f1ab58SBarry Smith 371b2f1ab58SBarry Smith /* 372b2f1ab58SBarry Smith spbas_delete : de-allocate the arrays owned by this matrix 373b2f1ab58SBarry Smith */ 374b2f1ab58SBarry Smith #undef __FUNCT__ 375b2f1ab58SBarry Smith #define __FUNCT__ "spbas_delete" 376b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix matrix) 377b2f1ab58SBarry Smith { 378b2f1ab58SBarry Smith PetscInt i; 379b2f1ab58SBarry Smith PetscErrorCode ierr; 3809ccc4a9bSBarry Smith 381b2f1ab58SBarry Smith PetscFunctionBegin; 3829ccc4a9bSBarry Smith if (matrix.block_data) { 383cb9801acSJed Brown ierr=PetscFree(matrix.alloc_icol);CHKERRQ(ierr); 384b2f1ab58SBarry Smith if (matrix.values) {ierr=PetscFree(matrix.alloc_val);CHKERRQ(ierr);} 3859ccc4a9bSBarry Smith } else { 3869ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.icols[i]);CHKERRQ(ierr);} 387b2f1ab58SBarry Smith ierr = PetscFree(matrix.icols);CHKERRQ(ierr); 3889ccc4a9bSBarry Smith if (matrix.values) { 3899ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.values[i]);CHKERRQ(ierr);} 390b2f1ab58SBarry Smith } 391b2f1ab58SBarry Smith } 392b2f1ab58SBarry Smith 393b2f1ab58SBarry Smith ierr=PetscFree(matrix.row_nnz);CHKERRQ(ierr); 394b2f1ab58SBarry Smith ierr=PetscFree(matrix.icols);CHKERRQ(ierr); 3959ccc4a9bSBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) {ierr=PetscFree(matrix.icol0);CHKERRQ(ierr);} 396c31cb41cSBarry Smith ierr=PetscFree(matrix.values);CHKERRQ(ierr); 397b2f1ab58SBarry Smith PetscFunctionReturn(0); 398b2f1ab58SBarry Smith } 399b2f1ab58SBarry Smith 400b2f1ab58SBarry Smith /* 401b2f1ab58SBarry Smith spbas_matrix_to_crs: 402b2f1ab58SBarry Smith Convert an spbas_matrix to compessed row storage 403b2f1ab58SBarry Smith */ 404b2f1ab58SBarry Smith #undef __FUNCT__ 405b2f1ab58SBarry Smith #define __FUNCT__ "spbas_matrix_to_crs" 406b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out) 407b2f1ab58SBarry Smith { 408b2f1ab58SBarry Smith PetscInt nrows = matrix_A.nrows; 409b2f1ab58SBarry Smith PetscInt nnz = matrix_A.nnz; 410b2f1ab58SBarry Smith PetscInt i,j,r_nnz,i0; 411b2f1ab58SBarry Smith PetscInt *irow; 412b2f1ab58SBarry Smith PetscInt *icol; 413b2f1ab58SBarry Smith PetscInt *icol_A; 414b2f1ab58SBarry Smith MatScalar *val; 415b2f1ab58SBarry Smith PetscScalar *val_A; 416b2f1ab58SBarry Smith PetscInt col_idx_type = matrix_A.col_idx_type; 417ace3abfcSBarry Smith PetscBool do_values = matrix_A.values ? PETSC_TRUE : PETSC_FALSE; 418b2f1ab58SBarry Smith PetscErrorCode ierr; 419b2f1ab58SBarry Smith 420b2f1ab58SBarry Smith PetscFunctionBegin; 421b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscInt) * (nrows+1), &irow);CHKERRQ(ierr); 422b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscInt) * nnz, &icol);CHKERRQ(ierr); 4239ccc4a9bSBarry Smith *icol_out = icol; 4249ccc4a9bSBarry Smith *irow_out = irow; 4259ccc4a9bSBarry Smith if (do_values) { 426b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(MatScalar) * nnz, &val);CHKERRQ(ierr); 427b2f1ab58SBarry Smith *val_out = val; *icol_out = icol; *irow_out=irow; 428b2f1ab58SBarry Smith } 429b2f1ab58SBarry Smith 430b2f1ab58SBarry Smith irow[0]=0; 4319ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 432b2f1ab58SBarry Smith r_nnz = matrix_A.row_nnz[i]; 433b2f1ab58SBarry Smith i0 = irow[i]; 434b2f1ab58SBarry Smith irow[i+1] = i0 + r_nnz; 435b2f1ab58SBarry Smith icol_A = matrix_A.icols[i]; 436b2f1ab58SBarry Smith 4379ccc4a9bSBarry Smith if (do_values) { 438b2f1ab58SBarry Smith val_A = matrix_A.values[i]; 4399ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 440b2f1ab58SBarry Smith icol[i0+j] = icol_A[j]; 441b2f1ab58SBarry Smith val[i0+j] = val_A[j]; 442b2f1ab58SBarry Smith } 4439ccc4a9bSBarry Smith } else { 444*2205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] = icol_A[j]; 445b2f1ab58SBarry Smith } 446b2f1ab58SBarry Smith 4479ccc4a9bSBarry Smith if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 448*2205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] += i; 449*2205254eSKarl Rupp } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 450b2f1ab58SBarry Smith i0 = matrix_A.icol0[i]; 451*2205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] += i0; 452b2f1ab58SBarry Smith } 453b2f1ab58SBarry Smith } 454b2f1ab58SBarry Smith PetscFunctionReturn(0); 455b2f1ab58SBarry Smith } 456b2f1ab58SBarry Smith 457b2f1ab58SBarry Smith 458b2f1ab58SBarry Smith /* 459b2f1ab58SBarry Smith spbas_transpose 460b2f1ab58SBarry Smith return the transpose of a matrix 461b2f1ab58SBarry Smith */ 462b2f1ab58SBarry Smith #undef __FUNCT__ 463b2f1ab58SBarry Smith #define __FUNCT__ "spbas_transpose" 464b2f1ab58SBarry Smith PetscErrorCode spbas_transpose(spbas_matrix in_matrix, spbas_matrix * result) 465b2f1ab58SBarry Smith { 466b2f1ab58SBarry Smith PetscInt col_idx_type = in_matrix.col_idx_type; 467b2f1ab58SBarry Smith PetscInt nnz = in_matrix.nnz; 468b2f1ab58SBarry Smith PetscInt ncols = in_matrix.nrows; 469b2f1ab58SBarry Smith PetscInt nrows = in_matrix.ncols; 470b2f1ab58SBarry Smith PetscInt i,j,k; 471b2f1ab58SBarry Smith PetscInt r_nnz; 472b2f1ab58SBarry Smith PetscInt *irow; 4734efc9174SBarry Smith PetscInt icol0 = 0; 474b2f1ab58SBarry Smith PetscScalar * val; 475b2f1ab58SBarry Smith PetscErrorCode ierr; 4764e1ff37aSBarry Smith 477b2f1ab58SBarry Smith PetscFunctionBegin; 478c328eeadSBarry Smith /* Copy input values */ 479b2f1ab58SBarry Smith result->nrows = nrows; 480b2f1ab58SBarry Smith result->ncols = ncols; 481b2f1ab58SBarry Smith result->nnz = nnz; 482b2f1ab58SBarry Smith result->col_idx_type = SPBAS_COLUMN_NUMBERS; 483b2f1ab58SBarry Smith result->block_data = PETSC_TRUE; 484b2f1ab58SBarry Smith 485c328eeadSBarry Smith /* Allocate sparseness pattern */ 48671d55d18SBarry Smith ierr = spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 487b2f1ab58SBarry Smith 488c328eeadSBarry Smith /* Count the number of nonzeros in each row */ 489*2205254eSKarl Rupp for (i = 0; i<nrows; i++) result->row_nnz[i] = 0; 490b2f1ab58SBarry Smith 4919ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 492b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 493b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 4949ccc4a9bSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 495*2205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[irow[j]]++; 4969ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 497*2205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[i+irow[j]]++; 4989ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 499b2f1ab58SBarry Smith icol0=in_matrix.icol0[i]; 500*2205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[icol0+irow[j]]++; 501b2f1ab58SBarry Smith } 502b2f1ab58SBarry Smith } 503b2f1ab58SBarry Smith 504c328eeadSBarry Smith /* Set the pointers to the data */ 505b2f1ab58SBarry Smith ierr = spbas_allocate_data(result);CHKERRQ(ierr); 506b2f1ab58SBarry Smith 507c328eeadSBarry Smith /* Reset the number of nonzeros in each row */ 508*2205254eSKarl Rupp for (i = 0; i<nrows; i++) result->row_nnz[i] = 0; 509b2f1ab58SBarry Smith 510c328eeadSBarry Smith /* Fill the data arrays */ 5119ccc4a9bSBarry Smith if (in_matrix.values) { 5129ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 513b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 514b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 515b2f1ab58SBarry Smith val = in_matrix.values[i]; 516b2f1ab58SBarry Smith 517*2205254eSKarl Rupp if (col_idx_type == SPBAS_COLUMN_NUMBERS) icol0 = 0; 518*2205254eSKarl Rupp else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0 = i; 519*2205254eSKarl Rupp else if (col_idx_type == SPBAS_OFFSET_ARRAY) icol0 = in_matrix.icol0[i]; 5209ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 521b2f1ab58SBarry Smith k = icol0 + irow[j]; 522b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 523b2f1ab58SBarry Smith result->values[k][result->row_nnz[k]] = val[j]; 524b2f1ab58SBarry Smith result->row_nnz[k]++; 525b2f1ab58SBarry Smith } 526b2f1ab58SBarry Smith } 5279ccc4a9bSBarry Smith } else { 5289ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 529b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 530b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 531b2f1ab58SBarry Smith 532*2205254eSKarl Rupp if (col_idx_type == SPBAS_COLUMN_NUMBERS) icol0=0; 533*2205254eSKarl Rupp else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0=i; 534*2205254eSKarl Rupp else if (col_idx_type == SPBAS_OFFSET_ARRAY) icol0=in_matrix.icol0[i]; 535b2f1ab58SBarry Smith 5369ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 537b2f1ab58SBarry Smith k = icol0 + irow[j]; 538b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 539b2f1ab58SBarry Smith result->row_nnz[k]++; 540b2f1ab58SBarry Smith } 541b2f1ab58SBarry Smith } 542b2f1ab58SBarry Smith } 543b2f1ab58SBarry Smith PetscFunctionReturn(0); 544b2f1ab58SBarry Smith } 545b2f1ab58SBarry Smith 546b2f1ab58SBarry Smith /* 547b2f1ab58SBarry Smith spbas_mergesort 548b2f1ab58SBarry Smith 549b2f1ab58SBarry Smith mergesort for an array of intergers and an array of associated 550b2f1ab58SBarry Smith reals 551b2f1ab58SBarry Smith 552b2f1ab58SBarry Smith on output, icol[0..nnz-1] is increasing; 553b2f1ab58SBarry Smith val[0..nnz-1] has undergone the same permutation as icol 554b2f1ab58SBarry Smith 555c328eeadSBarry Smith NB: val may be PETSC_NULL: in that case, only the integers are sorted 556b2f1ab58SBarry Smith 557b2f1ab58SBarry Smith */ 558b2f1ab58SBarry Smith #undef __FUNCT__ 559b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort" 560b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val) 561b2f1ab58SBarry Smith { 562c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 563c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 564c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 565c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 566c328eeadSBarry Smith PetscScalar *valloc=PETSC_NULL; 567c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 568b2f1ab58SBarry Smith PetscScalar *vswap; 569c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 570c328eeadSBarry Smith PetscScalar *vhlp1=PETSC_NULL; /* (arrays under construction) */ 571c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 572c328eeadSBarry Smith PetscScalar *vhlp2=PETSC_NULL; 573b2f1ab58SBarry Smith PetscErrorCode ierr; 574b2f1ab58SBarry Smith 575b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&ialloc);CHKERRQ(ierr); 576b2f1ab58SBarry Smith ihlp1 = ialloc; 577b2f1ab58SBarry Smith ihlp2 = icol; 578b2f1ab58SBarry Smith 5799ccc4a9bSBarry Smith if (val) { 580b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscScalar),&valloc);CHKERRQ(ierr); 581b2f1ab58SBarry Smith vhlp1 = valloc; 582b2f1ab58SBarry Smith vhlp2 = val; 583b2f1ab58SBarry Smith } 584b2f1ab58SBarry Smith 585b2f1ab58SBarry Smith 586c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 5879ccc4a9bSBarry Smith for (istep=1; istep<nnz; istep*=2) { 588c328eeadSBarry Smith /* 589c328eeadSBarry Smith Combine sorted parts 590c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 591c328eeadSBarry Smith of ihlp2 and vhlp2 592c328eeadSBarry Smith 593c328eeadSBarry Smith into one sorted part 594c328eeadSBarry Smith istart:istart+2*istep-1 595c328eeadSBarry Smith of ihlp1 and vhlp1 596c328eeadSBarry Smith */ 5979ccc4a9bSBarry Smith for (istart=0; istart<nnz; istart+=2*istep) { 598c328eeadSBarry Smith /* Set counters and bound array part endings */ 599*2205254eSKarl Rupp i1=istart; i1end = i1+istep; if (i1end>nnz) i1end=nnz; 600*2205254eSKarl Rupp i2=istart+istep; i2end = i2+istep; if (i2end>nnz) i2end=nnz; 601b2f1ab58SBarry Smith 602c328eeadSBarry Smith /* Merge the two array parts */ 6039ccc4a9bSBarry Smith if (val) { 6049ccc4a9bSBarry Smith for (i=istart; i<i2end; i++) { 6059ccc4a9bSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 606b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 607b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 608b2f1ab58SBarry Smith i1++; 6099ccc4a9bSBarry Smith } else if (i2<i2end) { 610b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 611b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i2]; 612b2f1ab58SBarry Smith i2++; 6139ccc4a9bSBarry Smith } else { 614b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 615b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 616b2f1ab58SBarry Smith i1++; 617b2f1ab58SBarry Smith } 618b2f1ab58SBarry Smith } 6199ccc4a9bSBarry Smith } else { 6206322f4bdSBarry Smith for (i=istart; i<i2end; i++) { 6216322f4bdSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 622b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 623b2f1ab58SBarry Smith i1++; 6246322f4bdSBarry Smith } else if (i2<i2end) { 625b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 626b2f1ab58SBarry Smith i2++; 6276322f4bdSBarry Smith } else { 628b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 629b2f1ab58SBarry Smith i1++; 630b2f1ab58SBarry Smith } 631b2f1ab58SBarry Smith } 632b2f1ab58SBarry Smith } 633b2f1ab58SBarry Smith } 634b2f1ab58SBarry Smith 635c328eeadSBarry Smith /* Swap the two array sets */ 636b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 637b2f1ab58SBarry Smith vswap = vhlp2; vhlp2 = vhlp1; vhlp1 = vswap; 638b2f1ab58SBarry Smith } 639b2f1ab58SBarry Smith 640c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 6416322f4bdSBarry Smith if (ihlp2 != icol) { 642*2205254eSKarl Rupp for (i=0; i<nnz; i++) icol[i] = ihlp2[i]; 6436322f4bdSBarry Smith if (val) { 644*2205254eSKarl Rupp for (i=0; i<nnz; i++) val[i] = vhlp2[i]; 645b2f1ab58SBarry Smith } 646b2f1ab58SBarry Smith } 647b2f1ab58SBarry Smith 648b2f1ab58SBarry Smith ierr = PetscFree(ialloc);CHKERRQ(ierr); 649b2f1ab58SBarry Smith if (val) {ierr = PetscFree(valloc);CHKERRQ(ierr);} 650b2f1ab58SBarry Smith PetscFunctionReturn(0); 651b2f1ab58SBarry Smith } 652b2f1ab58SBarry Smith 653b2f1ab58SBarry Smith /* 654b2f1ab58SBarry Smith spbas_apply_reordering_rows: 655b2f1ab58SBarry Smith apply the given reordering to the rows: matrix_A = matrix_A(perm,:); 656b2f1ab58SBarry Smith */ 657b2f1ab58SBarry Smith #undef __FUNCT__ 658b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_rows" 659b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_rows(spbas_matrix *matrix_A, const PetscInt *permutation) 660b2f1ab58SBarry Smith { 661b2f1ab58SBarry Smith PetscInt i,j,ip; 662b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 663b2f1ab58SBarry Smith PetscInt * row_nnz; 664b2f1ab58SBarry Smith PetscInt **icols; 665ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 666c328eeadSBarry Smith PetscScalar **vals = PETSC_NULL; 667b2f1ab58SBarry Smith PetscErrorCode ierr; 668b2f1ab58SBarry Smith 669b2f1ab58SBarry Smith PetscFunctionBegin; 670e32f2f54SBarry 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"); 671b2f1ab58SBarry Smith 6726322f4bdSBarry Smith if (do_values) { 673b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscScalar*)*nrows, &vals);CHKERRQ(ierr); 674b2f1ab58SBarry Smith } 675b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*nrows, &row_nnz);CHKERRQ(ierr); 676b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscInt*)*nrows, &icols);CHKERRQ(ierr); 677b2f1ab58SBarry Smith 6786322f4bdSBarry Smith for (i=0; i<nrows;i++) { 679b2f1ab58SBarry Smith ip = permutation[i]; 680*2205254eSKarl Rupp if (do_values) vals[i] = matrix_A->values[ip]; 681b2f1ab58SBarry Smith icols[i] = matrix_A->icols[ip]; 682b2f1ab58SBarry Smith row_nnz[i] = matrix_A->row_nnz[ip]; 683*2205254eSKarl Rupp for (j=0; j<row_nnz[i]; j++) icols[i][j] += ip-i; 684b2f1ab58SBarry Smith } 685b2f1ab58SBarry Smith 686b2f1ab58SBarry Smith if (do_values) { ierr = PetscFree(matrix_A->values);CHKERRQ(ierr);} 687b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->icols);CHKERRQ(ierr); 688b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->row_nnz);CHKERRQ(ierr); 689b2f1ab58SBarry Smith 690*2205254eSKarl Rupp if (do_values) matrix_A->values = vals; 691b2f1ab58SBarry Smith matrix_A->icols = icols; 692b2f1ab58SBarry Smith matrix_A->row_nnz = row_nnz; 693b2f1ab58SBarry Smith PetscFunctionReturn(0); 694b2f1ab58SBarry Smith } 695b2f1ab58SBarry Smith 696b2f1ab58SBarry Smith 697b2f1ab58SBarry Smith /* 698b2f1ab58SBarry Smith spbas_apply_reordering_cols: 699b2f1ab58SBarry Smith apply the given reordering to the columns: matrix_A(:,perm) = matrix_A; 700b2f1ab58SBarry Smith */ 701b2f1ab58SBarry Smith #undef __FUNCT__ 702b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_cols" 703b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_cols(spbas_matrix *matrix_A,const PetscInt *permutation) 704b2f1ab58SBarry Smith { 705b2f1ab58SBarry Smith PetscInt i,j; 706b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 707b2f1ab58SBarry Smith PetscInt row_nnz; 708b2f1ab58SBarry Smith PetscInt *icols; 709ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 710c328eeadSBarry Smith PetscScalar *vals = PETSC_NULL; 711b2f1ab58SBarry Smith PetscErrorCode ierr; 712b2f1ab58SBarry Smith 713b2f1ab58SBarry Smith PetscFunctionBegin; 714e32f2f54SBarry 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"); 715b2f1ab58SBarry Smith 7166322f4bdSBarry Smith for (i=0; i<nrows;i++) { 717b2f1ab58SBarry Smith icols = matrix_A->icols[i]; 718b2f1ab58SBarry Smith row_nnz = matrix_A->row_nnz[i]; 719*2205254eSKarl Rupp if (do_values) vals = matrix_A->values[i]; 720b2f1ab58SBarry Smith 7216322f4bdSBarry Smith for (j=0; j<row_nnz; j++) { 722b2f1ab58SBarry Smith icols[j] = permutation[i+icols[j]]-i; 723b2f1ab58SBarry Smith } 724b2f1ab58SBarry Smith ierr = spbas_mergesort(row_nnz, icols, vals);CHKERRQ(ierr); 725b2f1ab58SBarry Smith } 726b2f1ab58SBarry Smith PetscFunctionReturn(0); 727b2f1ab58SBarry Smith } 728b2f1ab58SBarry Smith 729b2f1ab58SBarry Smith /* 730b2f1ab58SBarry Smith spbas_apply_reordering: 731b2f1ab58SBarry Smith apply the given reordering: matrix_A(perm,perm) = matrix_A; 732b2f1ab58SBarry Smith */ 733b2f1ab58SBarry Smith #undef __FUNCT__ 734b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering" 735b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering(spbas_matrix *matrix_A, const PetscInt *permutation, const PetscInt * inv_perm) 736b2f1ab58SBarry Smith { 737b2f1ab58SBarry Smith PetscErrorCode ierr; 7385fd66863SKarl Rupp 739b2f1ab58SBarry Smith PetscFunctionBegin; 740b2f1ab58SBarry Smith ierr = spbas_apply_reordering_rows(matrix_A, inv_perm);CHKERRQ(ierr); 741b2f1ab58SBarry Smith ierr = spbas_apply_reordering_cols(matrix_A, permutation);CHKERRQ(ierr); 742b2f1ab58SBarry Smith PetscFunctionReturn(0); 743b2f1ab58SBarry Smith } 744b2f1ab58SBarry Smith 745b2f1ab58SBarry Smith #undef __FUNCT__ 746b2f1ab58SBarry Smith #define __FUNCT__ "spbas_pattern_only" 747b2f1ab58SBarry Smith PetscErrorCode spbas_pattern_only(PetscInt nrows, PetscInt ncols, PetscInt *ai, PetscInt *aj, spbas_matrix * result) 748b2f1ab58SBarry Smith { 749b2f1ab58SBarry Smith spbas_matrix retval; 750b2f1ab58SBarry Smith PetscInt i, j, i0, r_nnz; 751b2f1ab58SBarry Smith PetscErrorCode ierr; 752b2f1ab58SBarry Smith 753b2f1ab58SBarry Smith PetscFunctionBegin; 754c328eeadSBarry Smith /* Copy input values */ 755b2f1ab58SBarry Smith retval.nrows = nrows; 756b2f1ab58SBarry Smith retval.ncols = ncols; 757b2f1ab58SBarry Smith retval.nnz = ai[nrows]; 758b2f1ab58SBarry Smith 759b2f1ab58SBarry Smith retval.block_data = PETSC_TRUE; 760b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 761b2f1ab58SBarry Smith 762c328eeadSBarry Smith /* Allocate output matrix */ 763b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(&retval, PETSC_FALSE);CHKERRQ(ierr); 764*2205254eSKarl Rupp for (i=0; i<nrows; i++) retval.row_nnz[i] = ai[i+1]-ai[i]; 765b2f1ab58SBarry Smith ierr = spbas_allocate_data(&retval);CHKERRQ(ierr); 766c328eeadSBarry Smith /* Copy the structure */ 7676322f4bdSBarry Smith for (i = 0; i<retval.nrows; i++) { 768b2f1ab58SBarry Smith i0 = ai[i]; 769b2f1ab58SBarry Smith r_nnz = ai[i+1]-i0; 770b2f1ab58SBarry Smith 7716322f4bdSBarry Smith for (j=0; j<r_nnz; j++) { 772b2f1ab58SBarry Smith retval.icols[i][j] = aj[i0+j]-i; 773b2f1ab58SBarry Smith } 774b2f1ab58SBarry Smith } 775b2f1ab58SBarry Smith *result = retval; 776b2f1ab58SBarry Smith PetscFunctionReturn(0); 777b2f1ab58SBarry Smith } 778b2f1ab58SBarry Smith 779b2f1ab58SBarry Smith 780b2f1ab58SBarry Smith /* 781b2f1ab58SBarry Smith spbas_mark_row_power: 782b2f1ab58SBarry Smith Mark the columns in row 'row' which are nonzero in 783b2f1ab58SBarry Smith matrix^2log(marker). 784b2f1ab58SBarry Smith */ 785b2f1ab58SBarry Smith #undef __FUNCT__ 786b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mark_row_power" 787be332245SKarl Rupp PetscErrorCode spbas_mark_row_power(PetscInt *iwork, /* marker-vector */ 788c328eeadSBarry Smith PetscInt row, /* row for which the columns are marked */ 789c328eeadSBarry Smith spbas_matrix * in_matrix, /* matrix for which the power is being calculated */ 790c328eeadSBarry Smith PetscInt marker, /* marker-value: 2^power */ 791c328eeadSBarry Smith PetscInt minmrk, /* lower bound for marked points */ 792c328eeadSBarry Smith PetscInt maxmrk) /* upper bound for marked points */ 793b2f1ab58SBarry Smith { 794b2f1ab58SBarry Smith PetscErrorCode ierr; 795b2f1ab58SBarry Smith PetscInt i,j, nnz; 796b2f1ab58SBarry Smith 797b2f1ab58SBarry Smith PetscFunctionBegin; 798b2f1ab58SBarry Smith nnz = in_matrix->row_nnz[row]; 799b2f1ab58SBarry Smith 800c328eeadSBarry Smith /* For higher powers, call this function recursively */ 8016322f4bdSBarry Smith if (marker>1) { 8026322f4bdSBarry Smith for (i=0; i<nnz; i++) { 803b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 8046322f4bdSBarry Smith if (minmrk<=j && j<maxmrk && iwork[j] < marker) { 80571d55d18SBarry Smith ierr = spbas_mark_row_power(iwork, row + in_matrix->icols[row][i],in_matrix, marker/2,minmrk,maxmrk);CHKERRQ(ierr); 806b2f1ab58SBarry Smith iwork[j] |= marker; 807b2f1ab58SBarry Smith } 808b2f1ab58SBarry Smith } 8096322f4bdSBarry Smith } else { 810c328eeadSBarry Smith /* Mark the columns reached */ 8116322f4bdSBarry Smith for (i=0; i<nnz; i++) { 812b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 813*2205254eSKarl Rupp if (minmrk<=j && j<maxmrk) iwork[j] |= 1; 814b2f1ab58SBarry Smith } 815b2f1ab58SBarry Smith } 816b2f1ab58SBarry Smith PetscFunctionReturn(0); 817b2f1ab58SBarry Smith } 818b2f1ab58SBarry Smith 819b2f1ab58SBarry Smith 820b2f1ab58SBarry Smith /* 821b2f1ab58SBarry Smith spbas_power 822b2f1ab58SBarry Smith Calculate sparseness patterns for incomplete Cholesky decompositions 823b2f1ab58SBarry Smith of a given order: (almost) all nonzeros of the matrix^(order+1) which 824b2f1ab58SBarry Smith are inside the band width are found and stored in the output sparseness 825b2f1ab58SBarry Smith pattern. 826b2f1ab58SBarry Smith */ 827b2f1ab58SBarry Smith #undef __FUNCT__ 828b2f1ab58SBarry Smith #define __FUNCT__ "spbas_power" 829b2f1ab58SBarry Smith PetscErrorCode spbas_power(spbas_matrix in_matrix,PetscInt power, spbas_matrix * result) 830b2f1ab58SBarry Smith { 831b2f1ab58SBarry Smith spbas_matrix retval; 832b2f1ab58SBarry Smith PetscInt nrows = in_matrix.nrows; 833b2f1ab58SBarry Smith PetscInt ncols = in_matrix.ncols; 834b2f1ab58SBarry Smith PetscInt i, j, kend; 835b2f1ab58SBarry Smith PetscInt nnz, inz; 836b2f1ab58SBarry Smith PetscInt *iwork; 837b2f1ab58SBarry Smith PetscInt marker; 838b2f1ab58SBarry Smith PetscInt maxmrk=0; 839b2f1ab58SBarry Smith PetscErrorCode ierr; 840b2f1ab58SBarry Smith 841b2f1ab58SBarry Smith PetscFunctionBegin; 842e32f2f54SBarry 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"); 843e32f2f54SBarry Smith if (ncols != nrows) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Dimension error\n"); 844e32f2f54SBarry Smith if (in_matrix.values) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input array must be sparseness pattern (no values)"); 845e32f2f54SBarry Smith if (power<=0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Power must be 1 or up"); 846b2f1ab58SBarry Smith 847c328eeadSBarry Smith /* Copy input values*/ 848b2f1ab58SBarry Smith retval.nrows = ncols; 849b2f1ab58SBarry Smith retval.ncols = nrows; 850b2f1ab58SBarry Smith retval.nnz = 0; 851b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 852b2f1ab58SBarry Smith retval.block_data = PETSC_FALSE; 853b2f1ab58SBarry Smith 854c328eeadSBarry Smith /* Allocate sparseness pattern */ 85571d55d18SBarry Smith ierr = spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 856b2f1ab58SBarry Smith 857c328eeadSBarry Smith /* Allocate marker array */ 858b2f1ab58SBarry Smith ierr = PetscMalloc(nrows * sizeof(PetscInt), &iwork);CHKERRQ(ierr); 859b2f1ab58SBarry Smith 860c328eeadSBarry Smith /* Erase the pattern for this row */ 8614e1ff37aSBarry Smith ierr = PetscMemzero((void*) iwork, retval.nrows*sizeof(PetscInt));CHKERRQ(ierr); 862b2f1ab58SBarry Smith 863c328eeadSBarry Smith /* Calculate marker values */ 864*2205254eSKarl Rupp marker = 1; for (i=1; i<power; i++) marker*=2; 865b2f1ab58SBarry Smith 8666322f4bdSBarry Smith for (i=0; i<nrows; i++) { 867c328eeadSBarry Smith /* Calculate the pattern for each row */ 868b2f1ab58SBarry Smith 869b2f1ab58SBarry Smith nnz = in_matrix.row_nnz[i]; 870b2f1ab58SBarry Smith kend = i+in_matrix.icols[i][nnz-1]; 871*2205254eSKarl Rupp if (maxmrk<=kend) maxmrk=kend+1; 87271d55d18SBarry Smith ierr = spbas_mark_row_power(iwork, i, &in_matrix, marker, i, maxmrk);CHKERRQ(ierr); 873b2f1ab58SBarry Smith 874c328eeadSBarry Smith /* Count the columns*/ 875b2f1ab58SBarry Smith nnz = 0; 876*2205254eSKarl Rupp for (j=i; j<maxmrk; j++) nnz+= (iwork[j]!=0); 877b2f1ab58SBarry Smith 878c328eeadSBarry Smith /* Allocate the column indices */ 879b2f1ab58SBarry Smith retval.row_nnz[i] = nnz; 880b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&retval.icols[i]);CHKERRQ(ierr); 881b2f1ab58SBarry Smith 882c328eeadSBarry Smith /* Administrate the column indices */ 883b2f1ab58SBarry Smith inz = 0; 8846322f4bdSBarry Smith for (j=i; j<maxmrk; j++) { 8856322f4bdSBarry Smith if (iwork[j]) { 886b2f1ab58SBarry Smith retval.icols[i][inz] = j-i; 887b2f1ab58SBarry Smith inz++; 888b2f1ab58SBarry Smith iwork[j]=0; 889b2f1ab58SBarry Smith } 890b2f1ab58SBarry Smith } 891b2f1ab58SBarry Smith retval.nnz += nnz; 892b2f1ab58SBarry Smith }; 893b2f1ab58SBarry Smith ierr = PetscFree(iwork);CHKERRQ(ierr); 894b2f1ab58SBarry Smith *result = retval; 895b2f1ab58SBarry Smith PetscFunctionReturn(0); 896b2f1ab58SBarry Smith } 897b2f1ab58SBarry Smith 898b2f1ab58SBarry Smith 899b2f1ab58SBarry Smith 900b2f1ab58SBarry Smith /* 901b2f1ab58SBarry Smith spbas_keep_upper: 902b2f1ab58SBarry Smith remove the lower part of the matrix: keep the upper part 903b2f1ab58SBarry Smith */ 904b2f1ab58SBarry Smith #undef __FUNCT__ 905b2f1ab58SBarry Smith #define __FUNCT__ "spbas_keep_upper" 906b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper(spbas_matrix * inout_matrix) 907b2f1ab58SBarry Smith { 908b2f1ab58SBarry Smith PetscInt i, j; 909b2f1ab58SBarry Smith PetscInt jstart; 910b2f1ab58SBarry Smith 911b2f1ab58SBarry Smith PetscFunctionBegin; 912e32f2f54SBarry Smith if (inout_matrix->block_data) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Not yet for block data matrices\n"); 9136322f4bdSBarry Smith for (i=0; i<inout_matrix->nrows; i++) { 9146322f4bdSBarry Smith for (jstart=0; (jstart<inout_matrix->row_nnz[i]) && (inout_matrix->icols[i][jstart]<0); jstart++) {} 9156322f4bdSBarry Smith if (jstart>0) { 9166322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9176322f4bdSBarry Smith inout_matrix->icols[i][j] = inout_matrix->icols[i][j+jstart]; 918b2f1ab58SBarry Smith } 919b2f1ab58SBarry Smith 9206322f4bdSBarry Smith if (inout_matrix->values) { 9216322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9226322f4bdSBarry Smith inout_matrix->values[i][j] = inout_matrix->values[i][j+jstart]; 923b2f1ab58SBarry Smith } 924b2f1ab58SBarry Smith } 925b2f1ab58SBarry Smith 926b2f1ab58SBarry Smith inout_matrix->row_nnz[i] -= jstart; 927b2f1ab58SBarry Smith 9286322f4bdSBarry Smith inout_matrix->icols[i] = (PetscInt*) realloc((void*) inout_matrix->icols[i], inout_matrix->row_nnz[i]*sizeof(PetscInt)); 929b2f1ab58SBarry Smith 9306322f4bdSBarry Smith if (inout_matrix->values) { 9316322f4bdSBarry Smith inout_matrix->values[i] = (PetscScalar*) realloc((void*) inout_matrix->values[i], inout_matrix->row_nnz[i]*sizeof(PetscScalar)); 932b2f1ab58SBarry Smith } 933b2f1ab58SBarry Smith inout_matrix->nnz -= jstart; 934b2f1ab58SBarry Smith } 935b2f1ab58SBarry Smith } 936b2f1ab58SBarry Smith PetscFunctionReturn(0); 937b2f1ab58SBarry Smith } 938b2f1ab58SBarry Smith 939