1b2f1ab58SBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 2845006b9SBarry Smith #include "../src/mat/impls/aij/seq/bas/spbas.h" 3b2f1ab58SBarry Smith 4845006b9SBarry Smith /*MC 5845006b9SBarry Smith MAT_SOLVER_BAS - 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 { 29c328eeadSBarry Smith long int memreq = 6 * sizeof(PetscInt) + /* nrows, ncols, nnz, n_alloc_icol, 30c328eeadSBarry Smith n_alloc_val, col_idx_type */ 31c328eeadSBarry Smith sizeof(PetscTruth) + /* block_data */ 32c328eeadSBarry Smith sizeof(PetscScalar**) + /* values */ 33c328eeadSBarry Smith sizeof(PetscScalar*) + /* alloc_val */ 34c328eeadSBarry Smith 2 * sizeof(int**) + /* icols, icols0 */ 35c328eeadSBarry Smith 2 * sizeof(PetscInt*) + /* row_nnz, alloc_icol */ 36c328eeadSBarry Smith matrix.nrows * sizeof(PetscInt) + /* row_nnz[*] */ 37c328eeadSBarry Smith matrix.nrows * sizeof(PetscInt*); /* icols[*] */ 38b2f1ab58SBarry Smith 39c328eeadSBarry Smith /* icol0[*] */ 404e1ff37aSBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) { memreq += matrix.nrows * sizeof(PetscInt); } 41b2f1ab58SBarry Smith 42c328eeadSBarry Smith /* icols[*][*] */ 434e1ff37aSBarry Smith if (matrix.block_data) { memreq += matrix.n_alloc_icol * sizeof(PetscInt); } 444e1ff37aSBarry Smith else { memreq += matrix.nnz * sizeof(PetscInt); } 45b2f1ab58SBarry Smith 464e1ff37aSBarry Smith if (matrix.values) { 47c328eeadSBarry Smith memreq += matrix.nrows * sizeof(PetscScalar*); /* values[*] */ 48c328eeadSBarry Smith /* values[*][*] */ 494e1ff37aSBarry Smith if (matrix.block_data) { memreq += matrix.n_alloc_val * sizeof(PetscScalar); } 504e1ff37aSBarry Smith else { memreq += matrix.nnz * sizeof(PetscScalar); } 51b2f1ab58SBarry Smith } 52b2f1ab58SBarry Smith return memreq; 53b2f1ab58SBarry Smith } 54b2f1ab58SBarry Smith 55b2f1ab58SBarry Smith /* 56b2f1ab58SBarry Smith spbas_allocate_pattern: 57b2f1ab58SBarry Smith allocate the pattern arrays row_nnz, icols and optionally values 58b2f1ab58SBarry Smith */ 59b2f1ab58SBarry Smith #undef __FUNCT__ 60b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_pattern" 61b2f1ab58SBarry Smith PetscErrorCode spbas_allocate_pattern( spbas_matrix * result, PetscTruth do_values) 62b2f1ab58SBarry Smith { 63b2f1ab58SBarry Smith PetscErrorCode ierr; 64b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 65b2f1ab58SBarry Smith PetscInt col_idx_type = result->col_idx_type; 66b2f1ab58SBarry Smith 67b2f1ab58SBarry Smith PetscFunctionBegin; 68c328eeadSBarry Smith /* Allocate sparseness pattern */ 69b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->row_nnz);CHKERRQ(ierr); 70b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt*),&result->icols);CHKERRQ(ierr); 71b2f1ab58SBarry Smith 72c328eeadSBarry Smith /* If offsets are given wrt an array, create array */ 734e1ff37aSBarry Smith if (col_idx_type == SPBAS_OFFSET_ARRAY) { 74b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->icol0);CHKERRQ(ierr); 754e1ff37aSBarry Smith } else { 76c328eeadSBarry Smith result->icol0 = PETSC_NULL; 77b2f1ab58SBarry Smith } 78b2f1ab58SBarry Smith 79c328eeadSBarry Smith /* If values are given, allocate values array */ 804e1ff37aSBarry Smith if (do_values) { 81c328eeadSBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscScalar*),&result->values);CHKERRQ(ierr); 824e1ff37aSBarry Smith } else { 83c328eeadSBarry Smith result->values = PETSC_NULL; 84b2f1ab58SBarry Smith } 85b2f1ab58SBarry Smith PetscFunctionReturn(0); 86b2f1ab58SBarry Smith } 87b2f1ab58SBarry Smith 88b2f1ab58SBarry Smith /* 89b2f1ab58SBarry Smith spbas_allocate_data: 90b2f1ab58SBarry Smith in case of block_data: 91b2f1ab58SBarry Smith Allocate the data arrays alloc_icol and optionally alloc_val, 92b2f1ab58SBarry Smith set appropriate pointers from icols and values; 93b2f1ab58SBarry Smith in case of !block_data: 94b2f1ab58SBarry Smith Allocate the arrays icols[i] and optionally values[i] 95b2f1ab58SBarry Smith */ 96b2f1ab58SBarry Smith #undef __FUNCT__ 97b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_data" 98b2f1ab58SBarry Smith PetscErrorCode spbas_allocate_data( spbas_matrix * result) 99b2f1ab58SBarry Smith { 100b2f1ab58SBarry Smith PetscInt i; 101b2f1ab58SBarry Smith PetscInt nnz = result->nnz; 102b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 103b2f1ab58SBarry Smith PetscInt r_nnz; 104b2f1ab58SBarry Smith PetscErrorCode ierr; 105c328eeadSBarry Smith PetscTruth do_values = (result->values != PETSC_NULL) ? PETSC_TRUE : PETSC_FALSE; 106b2f1ab58SBarry Smith PetscTruth block_data = result->block_data; 107b2f1ab58SBarry Smith 108b2f1ab58SBarry Smith PetscFunctionBegin; 1094e1ff37aSBarry Smith if (block_data) { 110c328eeadSBarry Smith /* Allocate the column number array and point to it */ 111b2f1ab58SBarry Smith result->n_alloc_icol = nnz; 112c328eeadSBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt), &result->alloc_icol);CHKERRQ(ierr); 113b2f1ab58SBarry Smith result->icols[0]= result->alloc_icol; 1144e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 115b2f1ab58SBarry Smith result->icols[i] = result->icols[i-1] + result->row_nnz[i-1]; 116b2f1ab58SBarry Smith } 117b2f1ab58SBarry Smith 118c328eeadSBarry Smith /* Allocate the value array and point to it */ 1194e1ff37aSBarry Smith if (do_values) { 120b2f1ab58SBarry Smith result->n_alloc_val= nnz; 121c328eeadSBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscScalar), &result->alloc_val);CHKERRQ(ierr); 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]; 130c328eeadSBarry Smith ierr = PetscMalloc(r_nnz*sizeof(PetscInt), &result->icols[i]);CHKERRQ(ierr); 131b2f1ab58SBarry Smith } 1324e1ff37aSBarry Smith if (do_values) { 1334e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 134b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 1354e1ff37aSBarry Smith ierr = PetscMalloc(r_nnz*sizeof(PetscScalar), &result->values[i]);CHKERRQ(ierr); 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 #undef __FUNCT__ 150b2f1ab58SBarry Smith #define __FUNCT__ "spbas_row_order_icol" 151b2f1ab58SBarry Smith int spbas_row_order_icol(PetscInt i1, PetscInt i2, PetscInt *irow_in, PetscInt *icol_in,PetscInt col_idx_type) 152b2f1ab58SBarry Smith { 153b2f1ab58SBarry Smith PetscInt j; 154b2f1ab58SBarry Smith PetscInt nnz1 = irow_in[i1+1] - irow_in[i1]; 155b2f1ab58SBarry Smith PetscInt nnz2 = irow_in[i2+1] - irow_in[i2]; 156b2f1ab58SBarry Smith PetscInt * icol1 = &icol_in[irow_in[i1]]; 157b2f1ab58SBarry Smith PetscInt * icol2 = &icol_in[irow_in[i2]]; 158b2f1ab58SBarry Smith 159b2f1ab58SBarry Smith if (nnz1<nnz2) {return -1;} 160b2f1ab58SBarry Smith if (nnz1>nnz2) {return 1;} 161b2f1ab58SBarry Smith 1624e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 1634e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 164b2f1ab58SBarry Smith if (icol1[j]< icol2[j]) {return -1;} 165b2f1ab58SBarry Smith if (icol1[j]> icol2[j]) {return 1;} 166b2f1ab58SBarry Smith } 1674e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 1684e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 169b2f1ab58SBarry Smith if (icol1[j]-i1< icol2[j]-i2) {return -1;} 170b2f1ab58SBarry Smith if (icol1[j]-i1> icol2[j]-i2) {return 1;} 171b2f1ab58SBarry Smith } 1724e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 1734e1ff37aSBarry Smith for (j=1; j<nnz1; j++) { 174b2f1ab58SBarry Smith if (icol1[j]-icol1[0] < icol2[j]-icol2[0]) {return -1;} 175b2f1ab58SBarry Smith if (icol1[j]-icol1[0] > icol2[j]-icol2[0]) {return 1;} 176b2f1ab58SBarry Smith } 177b2f1ab58SBarry Smith } 178b2f1ab58SBarry Smith return 0; 179b2f1ab58SBarry Smith } 180b2f1ab58SBarry Smith 181b2f1ab58SBarry Smith /* 182b2f1ab58SBarry Smith spbas_mergesort_icols: 183b2f1ab58SBarry Smith return a sorting of the rows in which identical sparseness patterns are 184b2f1ab58SBarry Smith next to each other 185b2f1ab58SBarry Smith */ 186b2f1ab58SBarry Smith #undef __FUNCT__ 187b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort_icols" 188b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort_icols( PetscInt nrows, PetscInt * irow_in, PetscInt * icol_in,PetscInt col_idx_type, PetscInt *isort) 189b2f1ab58SBarry Smith { 190b2f1ab58SBarry Smith PetscErrorCode ierr; 191c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 192c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 193c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 194c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 195c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 196c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 197c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 198b2f1ab58SBarry Smith 199b2f1ab58SBarry Smith PetscFunctionBegin; 200b2f1ab58SBarry Smith 201b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&ialloc);CHKERRQ(ierr); 202b2f1ab58SBarry Smith 203b2f1ab58SBarry Smith ihlp1 = ialloc; 204b2f1ab58SBarry Smith ihlp2 = isort; 205b2f1ab58SBarry Smith 206c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 2074e1ff37aSBarry Smith for (istep=1; istep<nrows; istep*=2) { 208c328eeadSBarry Smith /* 209c328eeadSBarry Smith Combine sorted parts 210c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 211c328eeadSBarry Smith of ihlp2 and vhlp2 212c328eeadSBarry Smith 213c328eeadSBarry Smith into one sorted part 214c328eeadSBarry Smith istart:istart+2*istep-1 215c328eeadSBarry Smith of ihlp1 and vhlp1 216c328eeadSBarry Smith */ 2174e1ff37aSBarry Smith for (istart=0; istart<nrows; istart+=2*istep) { 218b2f1ab58SBarry Smith 219c328eeadSBarry Smith /* Set counters and bound array part endings */ 220b2f1ab58SBarry Smith i1=istart; i1end = i1+istep; if (i1end>nrows) {i1end=nrows;} 221b2f1ab58SBarry Smith i2=istart+istep; i2end = i2+istep; if (i2end>nrows) {i2end=nrows;} 222b2f1ab58SBarry Smith 223c328eeadSBarry Smith /* Merge the two array parts */ 2244e1ff37aSBarry Smith for (i=istart; i<i2end; i++) { 2254e1ff37aSBarry Smith if (i1<i1end && i2<i2end && spbas_row_order_icol( ihlp2[i1], ihlp2[i2], irow_in, icol_in, col_idx_type) < 0) { 226b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 227b2f1ab58SBarry Smith i1++; 2284e1ff37aSBarry Smith } else if (i2<i2end ) { 229b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 230b2f1ab58SBarry Smith i2++; 2314e1ff37aSBarry Smith } else { 232b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 233b2f1ab58SBarry Smith i1++; 234b2f1ab58SBarry Smith } 235b2f1ab58SBarry Smith } 236b2f1ab58SBarry Smith } 237b2f1ab58SBarry Smith 238c328eeadSBarry Smith /* Swap the two array sets */ 239b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 240b2f1ab58SBarry Smith } 241b2f1ab58SBarry Smith 242c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 2434e1ff37aSBarry Smith if (ihlp2 != isort) { 244b2f1ab58SBarry Smith for (i=0; i<nrows; i++) { isort[i] = ihlp2[i]; } 245b2f1ab58SBarry Smith } 246b2f1ab58SBarry Smith ierr = PetscFree(ialloc);CHKERRQ(ierr); 247b2f1ab58SBarry Smith PetscFunctionReturn(0); 248b2f1ab58SBarry Smith } 249b2f1ab58SBarry Smith 250b2f1ab58SBarry Smith 251b2f1ab58SBarry Smith 252b2f1ab58SBarry Smith /* 253b2f1ab58SBarry Smith spbas_compress_pattern: 254b2f1ab58SBarry Smith calculate a compressed sparseness pattern for a sparseness pattern 255b2f1ab58SBarry Smith given in compressed row storage. The compressed sparseness pattern may 256b2f1ab58SBarry Smith require (much) less memory. 257b2f1ab58SBarry Smith */ 258b2f1ab58SBarry Smith #undef __FUNCT__ 259b2f1ab58SBarry Smith #define __FUNCT__ "spbas_compress_pattern" 2604e1ff37aSBarry 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) 261b2f1ab58SBarry Smith { 262b2f1ab58SBarry Smith PetscInt nnz = irow_in[nrows]; 263b2f1ab58SBarry Smith long int mem_orig = (nrows + nnz) * sizeof(PetscInt); 264b2f1ab58SBarry Smith long int mem_compressed; 265b2f1ab58SBarry Smith PetscErrorCode ierr; 266b2f1ab58SBarry Smith PetscInt *isort; 267b2f1ab58SBarry Smith PetscInt *icols; 268b2f1ab58SBarry Smith PetscInt row_nnz; 269b2f1ab58SBarry Smith PetscInt *ipoint; 270b2f1ab58SBarry Smith PetscTruth *used; 271b2f1ab58SBarry Smith PetscInt ptr; 272b2f1ab58SBarry Smith PetscInt i,j; 273b2f1ab58SBarry Smith const PetscTruth no_values = PETSC_FALSE; 274b2f1ab58SBarry Smith 275b2f1ab58SBarry Smith PetscFunctionBegin; 276c328eeadSBarry Smith /* Allocate the structure of the new matrix */ 277b2f1ab58SBarry Smith B->nrows = nrows; 278b2f1ab58SBarry Smith B->ncols = ncols; 279b2f1ab58SBarry Smith B->nnz = nnz; 280b2f1ab58SBarry Smith B->col_idx_type= col_idx_type; 281b2f1ab58SBarry Smith B->block_data = PETSC_TRUE; 282b2f1ab58SBarry Smith ierr = spbas_allocate_pattern( B, no_values);CHKERRQ(ierr); 283b2f1ab58SBarry Smith 284c328eeadSBarry Smith /* When using an offset array, set it */ 2854e1ff37aSBarry Smith if (col_idx_type==SPBAS_OFFSET_ARRAY) { 286b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {B->icol0[i] = icol_in[irow_in[i]];} 287b2f1ab58SBarry Smith } 288b2f1ab58SBarry Smith 289c328eeadSBarry Smith /* Allocate the ordering for the rows */ 290b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&isort);CHKERRQ(ierr); 291b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&ipoint);CHKERRQ(ierr); 292b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscTruth),&used);CHKERRQ(ierr); 293b2f1ab58SBarry Smith 294c328eeadSBarry Smith /* Initialize the sorting */ 2954e1ff37aSBarry Smith ierr = PetscMemzero((void*) used, nrows*sizeof(PetscTruth));CHKERRQ(ierr); 2964e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 297b2f1ab58SBarry Smith B->row_nnz[i] = irow_in[i+1]-irow_in[i]; 298b2f1ab58SBarry Smith isort[i] = i; 299b2f1ab58SBarry Smith ipoint[i]= i; 300b2f1ab58SBarry Smith } 301b2f1ab58SBarry Smith 302c328eeadSBarry Smith /* Sort the rows so that identical columns will be next to each other */ 303b2f1ab58SBarry Smith ierr = spbas_mergesort_icols( nrows, irow_in, icol_in, col_idx_type, isort);CHKERRQ(ierr); 30471d55d18SBarry Smith ierr = PetscInfo(PETSC_NULL,"Rows have been sorted for patterns\n");CHKERRQ(ierr); 305b2f1ab58SBarry Smith 306c328eeadSBarry Smith /* Replace identical rows with the first one in the list */ 3074e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 3084e1ff37aSBarry Smith if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) { 309b2f1ab58SBarry Smith ipoint[isort[i]] = ipoint[isort[i-1]]; 310b2f1ab58SBarry Smith } 311b2f1ab58SBarry Smith } 312b2f1ab58SBarry Smith 313c328eeadSBarry Smith /* Collect the rows which are used*/ 314b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {used[ipoint[i]] = PETSC_TRUE;} 315b2f1ab58SBarry Smith 316c328eeadSBarry Smith /* Calculate needed memory */ 317b2f1ab58SBarry Smith B->n_alloc_icol = 0; 3184e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 319b2f1ab58SBarry Smith if (used[i]) {B->n_alloc_icol += B->row_nnz[i];} 320b2f1ab58SBarry Smith } 32171d55d18SBarry Smith ierr = PetscMalloc(B->n_alloc_icol*sizeof(PetscInt),&B->alloc_icol);CHKERRQ(ierr); 322b2f1ab58SBarry Smith 323c328eeadSBarry Smith /* Fill in the diagonal offsets for the rows which store their own data */ 324b2f1ab58SBarry Smith ptr = 0; 3254e1ff37aSBarry Smith for (i=0; i<B->nrows; i++) { 3264e1ff37aSBarry Smith if (used[i]) { 327b2f1ab58SBarry Smith B->icols[i] = &B->alloc_icol[ptr]; 328b2f1ab58SBarry Smith icols = &icol_in[irow_in[i]]; 329b2f1ab58SBarry Smith row_nnz = B->row_nnz[i]; 3304e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 3314e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 332b2f1ab58SBarry Smith B->icols[i][j] = icols[j]; 333b2f1ab58SBarry Smith } 3344e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 3354e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 336b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-i; 337b2f1ab58SBarry Smith } 3384e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 3394e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 340b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-icols[0]; 341b2f1ab58SBarry Smith } 342b2f1ab58SBarry Smith } 343b2f1ab58SBarry Smith ptr += B->row_nnz[i]; 344b2f1ab58SBarry Smith } 345b2f1ab58SBarry Smith } 346b2f1ab58SBarry Smith 347c328eeadSBarry Smith /* Point to the right places for all data */ 3484e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 349b2f1ab58SBarry Smith B->icols[i] = B->icols[ipoint[i]]; 350b2f1ab58SBarry Smith } 351c328eeadSBarry Smith ierr = PetscInfo(PETSC_NULL,"Row patterns have been compressed\n");CHKERRQ(ierr); 352c328eeadSBarry Smith ierr = PetscInfo1(PETSC_NULL," (%G nonzeros per row)\n", (PetscReal) nnz / (PetscReal) nrows);CHKERRQ(ierr); 353b2f1ab58SBarry Smith 354b2f1ab58SBarry Smith ierr=PetscFree(isort);CHKERRQ(ierr); 355b2f1ab58SBarry Smith ierr=PetscFree(used);CHKERRQ(ierr); 356b2f1ab58SBarry Smith ierr=PetscFree(ipoint);CHKERRQ(ierr); 357b2f1ab58SBarry Smith 358b2f1ab58SBarry Smith mem_compressed = spbas_memory_requirement( *B ); 3594e1ff37aSBarry Smith *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig; 360b2f1ab58SBarry Smith PetscFunctionReturn(0); 361b2f1ab58SBarry Smith } 362b2f1ab58SBarry Smith 363b2f1ab58SBarry Smith /* 364b2f1ab58SBarry Smith spbas_incomplete_cholesky 365b2f1ab58SBarry Smith Incomplete Cholesky decomposition 366b2f1ab58SBarry Smith */ 367845006b9SBarry Smith #include "../src/mat/impls/aij/seq/bas/spbas_cholesky.h" 368b2f1ab58SBarry Smith 369b2f1ab58SBarry Smith /* 370b2f1ab58SBarry Smith spbas_delete : de-allocate the arrays owned by this matrix 371b2f1ab58SBarry Smith */ 372b2f1ab58SBarry Smith #undef __FUNCT__ 373b2f1ab58SBarry Smith #define __FUNCT__ "spbas_delete" 374b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix matrix) 375b2f1ab58SBarry Smith { 376b2f1ab58SBarry Smith PetscInt i; 377b2f1ab58SBarry Smith PetscErrorCode ierr; 3789ccc4a9bSBarry Smith 379b2f1ab58SBarry Smith PetscFunctionBegin; 3809ccc4a9bSBarry Smith if (matrix.block_data) { 381b2f1ab58SBarry Smith ierr=PetscFree(matrix.alloc_icol);CHKERRQ(ierr) 382b2f1ab58SBarry Smith if (matrix.values){ierr=PetscFree(matrix.alloc_val);CHKERRQ(ierr);} 3839ccc4a9bSBarry Smith } else { 3849ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.icols[i]);CHKERRQ(ierr);} 385b2f1ab58SBarry Smith ierr = PetscFree(matrix.icols);CHKERRQ(ierr); 3869ccc4a9bSBarry Smith if (matrix.values) { 3879ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.values[i]);CHKERRQ(ierr);} 388b2f1ab58SBarry Smith } 389b2f1ab58SBarry Smith } 390b2f1ab58SBarry Smith 391b2f1ab58SBarry Smith ierr=PetscFree(matrix.row_nnz);CHKERRQ(ierr); 392b2f1ab58SBarry Smith ierr=PetscFree(matrix.icols);CHKERRQ(ierr); 3939ccc4a9bSBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) {ierr=PetscFree(matrix.icol0);CHKERRQ(ierr);} 3949ccc4a9bSBarry Smith if (matrix.values) {ierr=PetscFree(matrix.values);CHKERRQ(ierr);} 395b2f1ab58SBarry Smith PetscFunctionReturn(0); 396b2f1ab58SBarry Smith } 397b2f1ab58SBarry Smith 398b2f1ab58SBarry Smith /* 399b2f1ab58SBarry Smith spbas_matrix_to_crs: 400b2f1ab58SBarry Smith Convert an spbas_matrix to compessed row storage 401b2f1ab58SBarry Smith */ 402b2f1ab58SBarry Smith #undef __FUNCT__ 403b2f1ab58SBarry Smith #define __FUNCT__ "spbas_matrix_to_crs" 404b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out) 405b2f1ab58SBarry Smith { 406b2f1ab58SBarry Smith PetscInt nrows = matrix_A.nrows; 407b2f1ab58SBarry Smith PetscInt nnz = matrix_A.nnz; 408b2f1ab58SBarry Smith PetscInt i,j,r_nnz,i0; 409b2f1ab58SBarry Smith PetscInt *irow; 410b2f1ab58SBarry Smith PetscInt *icol; 411b2f1ab58SBarry Smith PetscInt *icol_A; 412b2f1ab58SBarry Smith MatScalar *val; 413b2f1ab58SBarry Smith PetscScalar *val_A; 414b2f1ab58SBarry Smith PetscInt col_idx_type = matrix_A.col_idx_type; 415d36aa34eSBarry Smith PetscTruth do_values = matrix_A.values ? PETSC_TRUE : PETSC_FALSE; 416b2f1ab58SBarry Smith PetscErrorCode ierr; 417b2f1ab58SBarry Smith 418b2f1ab58SBarry Smith PetscFunctionBegin; 419b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt) * (nrows+1), &irow);CHKERRQ(ierr); 420b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt) * nnz, &icol);CHKERRQ(ierr); 4219ccc4a9bSBarry Smith *icol_out = icol; 4229ccc4a9bSBarry Smith *irow_out=irow; 4239ccc4a9bSBarry Smith if (do_values) { 424b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(MatScalar) * nnz, &val);CHKERRQ(ierr); 425b2f1ab58SBarry Smith *val_out = val; *icol_out = icol; *irow_out=irow; 426b2f1ab58SBarry Smith } 427b2f1ab58SBarry Smith 428b2f1ab58SBarry Smith irow[0]=0; 4299ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 430b2f1ab58SBarry Smith r_nnz = matrix_A.row_nnz[i]; 431b2f1ab58SBarry Smith i0 = irow[i]; 432b2f1ab58SBarry Smith irow[i+1] = i0 + r_nnz; 433b2f1ab58SBarry Smith icol_A = matrix_A.icols[i]; 434b2f1ab58SBarry Smith 4359ccc4a9bSBarry Smith if (do_values) { 436b2f1ab58SBarry Smith val_A = matrix_A.values[i]; 4379ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 438b2f1ab58SBarry Smith icol[i0+j] = icol_A[j]; 439b2f1ab58SBarry Smith val[i0+j] = val_A[j]; 440b2f1ab58SBarry Smith } 4419ccc4a9bSBarry Smith } else { 442b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] = icol_A[j]; } 443b2f1ab58SBarry Smith } 444b2f1ab58SBarry Smith 4459ccc4a9bSBarry Smith if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 446b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] += i; } 447b2f1ab58SBarry Smith } 4489ccc4a9bSBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 449b2f1ab58SBarry Smith i0 = matrix_A.icol0[i]; 450b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] += i0; } 451b2f1ab58SBarry Smith } 452b2f1ab58SBarry Smith } 453b2f1ab58SBarry Smith PetscFunctionReturn(0); 454b2f1ab58SBarry Smith } 455b2f1ab58SBarry Smith 456b2f1ab58SBarry Smith 457b2f1ab58SBarry Smith /* 458b2f1ab58SBarry Smith spbas_transpose 459b2f1ab58SBarry Smith return the transpose of a matrix 460b2f1ab58SBarry Smith */ 461b2f1ab58SBarry Smith #undef __FUNCT__ 462b2f1ab58SBarry Smith #define __FUNCT__ "spbas_transpose" 463b2f1ab58SBarry Smith PetscErrorCode spbas_transpose( spbas_matrix in_matrix, spbas_matrix * result) 464b2f1ab58SBarry Smith { 465b2f1ab58SBarry Smith PetscInt col_idx_type = in_matrix.col_idx_type; 466b2f1ab58SBarry Smith PetscInt nnz = in_matrix.nnz; 467b2f1ab58SBarry Smith PetscInt ncols = in_matrix.nrows; 468b2f1ab58SBarry Smith PetscInt nrows = in_matrix.ncols; 469b2f1ab58SBarry Smith PetscInt i,j,k; 470b2f1ab58SBarry Smith PetscInt r_nnz; 471b2f1ab58SBarry Smith PetscInt *irow; 472*4efc9174SBarry Smith PetscInt icol0 = 0; 473b2f1ab58SBarry Smith PetscScalar * val; 474b2f1ab58SBarry Smith PetscErrorCode ierr; 4754e1ff37aSBarry Smith 476b2f1ab58SBarry Smith PetscFunctionBegin; 477c328eeadSBarry Smith /* Copy input values */ 478b2f1ab58SBarry Smith result->nrows = nrows; 479b2f1ab58SBarry Smith result->ncols = ncols; 480b2f1ab58SBarry Smith result->nnz = nnz; 481b2f1ab58SBarry Smith result->col_idx_type = SPBAS_COLUMN_NUMBERS; 482b2f1ab58SBarry Smith result->block_data = PETSC_TRUE; 483b2f1ab58SBarry Smith 484c328eeadSBarry Smith /* Allocate sparseness pattern */ 48571d55d18SBarry Smith ierr = spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 486b2f1ab58SBarry Smith 487c328eeadSBarry Smith /* Count the number of nonzeros in each row */ 488b2f1ab58SBarry Smith for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; } 489b2f1ab58SBarry Smith 4909ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 491b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 492b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 4939ccc4a9bSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 494b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[irow[j]]++; } 4959ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 496b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[i+irow[j]]++; } 4979ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 498b2f1ab58SBarry Smith icol0=in_matrix.icol0[i]; 499b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[icol0+irow[j]]++; } 500b2f1ab58SBarry Smith } 501b2f1ab58SBarry Smith } 502b2f1ab58SBarry Smith 503c328eeadSBarry Smith /* Set the pointers to the data */ 504b2f1ab58SBarry Smith ierr = spbas_allocate_data(result);CHKERRQ(ierr); 505b2f1ab58SBarry Smith 506c328eeadSBarry Smith /* Reset the number of nonzeros in each row */ 507b2f1ab58SBarry Smith for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; } 508b2f1ab58SBarry Smith 509c328eeadSBarry Smith /* Fill the data arrays */ 5109ccc4a9bSBarry Smith if (in_matrix.values) { 5119ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 512b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 513b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 514b2f1ab58SBarry Smith val = in_matrix.values[i]; 515b2f1ab58SBarry Smith 516b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) {icol0=0;} 517b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;} 518*4efc9174SBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) {icol0=in_matrix.icol0[i];} 5199ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 520b2f1ab58SBarry Smith k = icol0 + irow[j]; 521b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 522b2f1ab58SBarry Smith result->values[k][result->row_nnz[k]] = val[j]; 523b2f1ab58SBarry Smith result->row_nnz[k]++; 524b2f1ab58SBarry Smith } 525b2f1ab58SBarry Smith } 5269ccc4a9bSBarry Smith } else { 5279ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 528b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 529b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 530b2f1ab58SBarry Smith 531b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) {icol0=0;} 532b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;} 5339ccc4a9bSBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) {icol0=in_matrix.icol0[i];} 534b2f1ab58SBarry Smith 5359ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 536b2f1ab58SBarry Smith k = icol0 + irow[j]; 537b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 538b2f1ab58SBarry Smith result->row_nnz[k]++; 539b2f1ab58SBarry Smith } 540b2f1ab58SBarry Smith } 541b2f1ab58SBarry Smith } 542b2f1ab58SBarry Smith PetscFunctionReturn(0); 543b2f1ab58SBarry Smith } 544b2f1ab58SBarry Smith 545b2f1ab58SBarry Smith /* 546b2f1ab58SBarry Smith spbas_mergesort 547b2f1ab58SBarry Smith 548b2f1ab58SBarry Smith mergesort for an array of intergers and an array of associated 549b2f1ab58SBarry Smith reals 550b2f1ab58SBarry Smith 551b2f1ab58SBarry Smith on output, icol[0..nnz-1] is increasing; 552b2f1ab58SBarry Smith val[0..nnz-1] has undergone the same permutation as icol 553b2f1ab58SBarry Smith 554c328eeadSBarry Smith NB: val may be PETSC_NULL: in that case, only the integers are sorted 555b2f1ab58SBarry Smith 556b2f1ab58SBarry Smith */ 557b2f1ab58SBarry Smith #undef __FUNCT__ 558b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort" 559b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val) 560b2f1ab58SBarry Smith { 561c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 562c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 563c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 564c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 565c328eeadSBarry Smith PetscScalar *valloc=PETSC_NULL; 566c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 567b2f1ab58SBarry Smith PetscScalar *vswap; 568c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 569c328eeadSBarry Smith PetscScalar *vhlp1=PETSC_NULL; /* (arrays under construction) */ 570c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 571c328eeadSBarry Smith PetscScalar *vhlp2=PETSC_NULL; 572b2f1ab58SBarry Smith PetscErrorCode ierr; 573b2f1ab58SBarry Smith 574b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&ialloc);CHKERRQ(ierr); 575b2f1ab58SBarry Smith ihlp1 = ialloc; 576b2f1ab58SBarry Smith ihlp2 = icol; 577b2f1ab58SBarry Smith 5789ccc4a9bSBarry Smith if (val) { 579b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscScalar),&valloc);CHKERRQ(ierr); 580b2f1ab58SBarry Smith vhlp1 = valloc; 581b2f1ab58SBarry Smith vhlp2 = val; 582b2f1ab58SBarry Smith } 583b2f1ab58SBarry Smith 584b2f1ab58SBarry Smith 585c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 5869ccc4a9bSBarry Smith for (istep=1; istep<nnz; istep*=2) { 587c328eeadSBarry Smith /* 588c328eeadSBarry Smith Combine sorted parts 589c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 590c328eeadSBarry Smith of ihlp2 and vhlp2 591c328eeadSBarry Smith 592c328eeadSBarry Smith into one sorted part 593c328eeadSBarry Smith istart:istart+2*istep-1 594c328eeadSBarry Smith of ihlp1 and vhlp1 595c328eeadSBarry Smith */ 5969ccc4a9bSBarry Smith for (istart=0; istart<nnz; istart+=2*istep) { 597b2f1ab58SBarry Smith 598c328eeadSBarry Smith /* Set counters and bound array part endings */ 599b2f1ab58SBarry Smith i1=istart; i1end = i1+istep; if (i1end>nnz) {i1end=nnz;} 600b2f1ab58SBarry Smith 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) { 642b2f1ab58SBarry Smith for (i=0; i<nnz; i++) { icol[i] = ihlp2[i]; } 6436322f4bdSBarry Smith if (val) { 644b2f1ab58SBarry Smith 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; 665d36aa34eSBarry Smith PetscTruth do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 666c328eeadSBarry Smith PetscScalar **vals=PETSC_NULL; 667b2f1ab58SBarry Smith PetscErrorCode ierr; 668b2f1ab58SBarry Smith 669b2f1ab58SBarry Smith PetscFunctionBegin; 6706322f4bdSBarry Smith if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ( 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]; 680b2f1ab58SBarry Smith 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]; 683b2f1ab58SBarry Smith 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 690b2f1ab58SBarry Smith if (do_values) { matrix_A->values = vals; } 691b2f1ab58SBarry Smith matrix_A->icols = icols; 692b2f1ab58SBarry Smith matrix_A->row_nnz = row_nnz; 693b2f1ab58SBarry Smith 694b2f1ab58SBarry Smith PetscFunctionReturn(0); 695b2f1ab58SBarry Smith } 696b2f1ab58SBarry Smith 697b2f1ab58SBarry Smith 698b2f1ab58SBarry Smith /* 699b2f1ab58SBarry Smith spbas_apply_reordering_cols: 700b2f1ab58SBarry Smith apply the given reordering to the columns: matrix_A(:,perm) = matrix_A; 701b2f1ab58SBarry Smith */ 702b2f1ab58SBarry Smith #undef __FUNCT__ 703b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_cols" 704b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_cols( spbas_matrix *matrix_A,const PetscInt *permutation) 705b2f1ab58SBarry Smith { 706b2f1ab58SBarry Smith PetscInt i,j; 707b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 708b2f1ab58SBarry Smith PetscInt row_nnz; 709b2f1ab58SBarry Smith PetscInt *icols; 710d36aa34eSBarry Smith PetscTruth do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 711c328eeadSBarry Smith PetscScalar *vals=PETSC_NULL; 712b2f1ab58SBarry Smith PetscErrorCode ierr; 713b2f1ab58SBarry Smith 714b2f1ab58SBarry Smith PetscFunctionBegin; 7156322f4bdSBarry Smith if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ( PETSC_ERR_SUP_SYS, "must have diagonal offsets in pattern\n"); 716b2f1ab58SBarry Smith 7176322f4bdSBarry Smith for (i=0; i<nrows;i++) { 718b2f1ab58SBarry Smith icols = matrix_A->icols[i]; 719b2f1ab58SBarry Smith row_nnz = matrix_A->row_nnz[i]; 720b2f1ab58SBarry Smith if (do_values) { vals = matrix_A->values[i]; } 721b2f1ab58SBarry Smith 7226322f4bdSBarry Smith for (j=0; j<row_nnz; j++) { 723b2f1ab58SBarry Smith icols[j] = permutation[i+icols[j]]-i; 724b2f1ab58SBarry Smith } 725b2f1ab58SBarry Smith ierr = spbas_mergesort(row_nnz, icols, vals);CHKERRQ(ierr); 726b2f1ab58SBarry Smith } 727b2f1ab58SBarry Smith PetscFunctionReturn(0); 728b2f1ab58SBarry Smith } 729b2f1ab58SBarry Smith 730b2f1ab58SBarry Smith /* 731b2f1ab58SBarry Smith spbas_apply_reordering: 732b2f1ab58SBarry Smith apply the given reordering: matrix_A(perm,perm) = matrix_A; 733b2f1ab58SBarry Smith */ 734b2f1ab58SBarry Smith #undef __FUNCT__ 735b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering" 736b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering( spbas_matrix *matrix_A, const PetscInt *permutation, const PetscInt * inv_perm) 737b2f1ab58SBarry Smith { 738b2f1ab58SBarry Smith PetscErrorCode ierr; 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; 754b2f1ab58SBarry Smith 755c328eeadSBarry Smith /* Copy input values */ 756b2f1ab58SBarry Smith retval.nrows = nrows; 757b2f1ab58SBarry Smith retval.ncols = ncols; 758b2f1ab58SBarry Smith retval.nnz = ai[nrows]; 759b2f1ab58SBarry Smith 760b2f1ab58SBarry Smith retval.block_data = PETSC_TRUE; 761b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 762b2f1ab58SBarry Smith 763c328eeadSBarry Smith /* Allocate output matrix */ 764b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(&retval, PETSC_FALSE);CHKERRQ(ierr); 765b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {retval.row_nnz[i] = ai[i+1]-ai[i];} 766b2f1ab58SBarry Smith ierr = spbas_allocate_data(&retval);CHKERRQ(ierr); 767c328eeadSBarry Smith /* Copy the structure */ 7686322f4bdSBarry Smith for (i = 0; i<retval.nrows; i++) { 769b2f1ab58SBarry Smith i0 = ai[i]; 770b2f1ab58SBarry Smith r_nnz = ai[i+1]-i0; 771b2f1ab58SBarry Smith 7726322f4bdSBarry Smith for (j=0; j<r_nnz; j++) { 773b2f1ab58SBarry Smith retval.icols[i][j] = aj[i0+j]-i; 774b2f1ab58SBarry Smith } 775b2f1ab58SBarry Smith } 776b2f1ab58SBarry Smith *result = retval; 777b2f1ab58SBarry Smith PetscFunctionReturn(0); 778b2f1ab58SBarry Smith } 779b2f1ab58SBarry Smith 780b2f1ab58SBarry Smith 781b2f1ab58SBarry Smith /* 782b2f1ab58SBarry Smith spbas_mark_row_power: 783b2f1ab58SBarry Smith Mark the columns in row 'row' which are nonzero in 784b2f1ab58SBarry Smith matrix^2log(marker). 785b2f1ab58SBarry Smith */ 786b2f1ab58SBarry Smith #undef __FUNCT__ 787b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mark_row_power" 788b2f1ab58SBarry Smith PetscErrorCode spbas_mark_row_power( 789c328eeadSBarry Smith PetscInt *iwork, /* marker-vector */ 790c328eeadSBarry Smith PetscInt row, /* row for which the columns are marked */ 791c328eeadSBarry Smith spbas_matrix * in_matrix, /* matrix for which the power is being calculated */ 792c328eeadSBarry Smith PetscInt marker, /* marker-value: 2^power */ 793c328eeadSBarry Smith PetscInt minmrk, /* lower bound for marked points */ 794c328eeadSBarry Smith PetscInt maxmrk) /* upper bound for marked points */ 795b2f1ab58SBarry Smith { 796b2f1ab58SBarry Smith PetscErrorCode ierr; 797b2f1ab58SBarry Smith PetscInt i,j, nnz; 798b2f1ab58SBarry Smith 799b2f1ab58SBarry Smith PetscFunctionBegin; 800b2f1ab58SBarry Smith nnz = in_matrix->row_nnz[row]; 801b2f1ab58SBarry Smith 802c328eeadSBarry Smith /* For higher powers, call this function recursively */ 8036322f4bdSBarry Smith if (marker>1) { 8046322f4bdSBarry Smith for (i=0; i<nnz;i++) { 805b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 8066322f4bdSBarry Smith if (minmrk<=j && j<maxmrk && iwork[j] < marker ) { 80771d55d18SBarry Smith ierr = spbas_mark_row_power( iwork, row + in_matrix->icols[row][i],in_matrix, marker/2,minmrk,maxmrk);CHKERRQ(ierr); 808b2f1ab58SBarry Smith iwork[j] |= marker; 809b2f1ab58SBarry Smith } 810b2f1ab58SBarry Smith } 8116322f4bdSBarry Smith } else { 812c328eeadSBarry Smith /* Mark the columns reached */ 8136322f4bdSBarry Smith for (i=0; i<nnz;i++) { 814b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 815b2f1ab58SBarry Smith if (minmrk<=j && j<maxmrk) { iwork[j] |= 1; } 816b2f1ab58SBarry Smith } 817b2f1ab58SBarry Smith } 818b2f1ab58SBarry Smith PetscFunctionReturn(0); 819b2f1ab58SBarry Smith } 820b2f1ab58SBarry Smith 821b2f1ab58SBarry Smith 822b2f1ab58SBarry Smith /* 823b2f1ab58SBarry Smith spbas_power 824b2f1ab58SBarry Smith Calculate sparseness patterns for incomplete Cholesky decompositions 825b2f1ab58SBarry Smith of a given order: (almost) all nonzeros of the matrix^(order+1) which 826b2f1ab58SBarry Smith are inside the band width are found and stored in the output sparseness 827b2f1ab58SBarry Smith pattern. 828b2f1ab58SBarry Smith */ 829b2f1ab58SBarry Smith #undef __FUNCT__ 830b2f1ab58SBarry Smith #define __FUNCT__ "spbas_power" 831b2f1ab58SBarry Smith PetscErrorCode spbas_power (spbas_matrix in_matrix,PetscInt power, spbas_matrix * result) 832b2f1ab58SBarry Smith { 833b2f1ab58SBarry Smith spbas_matrix retval; 834b2f1ab58SBarry Smith PetscInt nrows = in_matrix.nrows; 835b2f1ab58SBarry Smith PetscInt ncols = in_matrix.ncols; 836b2f1ab58SBarry Smith PetscInt i, j, kend; 837b2f1ab58SBarry Smith PetscInt nnz, inz; 838b2f1ab58SBarry Smith PetscInt *iwork; 839b2f1ab58SBarry Smith PetscInt marker; 840b2f1ab58SBarry Smith PetscInt maxmrk=0; 841b2f1ab58SBarry Smith PetscErrorCode ierr; 842b2f1ab58SBarry Smith 843b2f1ab58SBarry Smith PetscFunctionBegin; 8446322f4bdSBarry Smith if (in_matrix.col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ( PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern\n"); 8456322f4bdSBarry Smith if (ncols != nrows) SETERRQ( PETSC_ERR_ARG_INCOMP, "Dimension error\n"); 8466322f4bdSBarry Smith if (in_matrix.values) SETERRQ( PETSC_ERR_ARG_INCOMP, "Input array must be sparseness pattern (no values)"); 8476322f4bdSBarry Smith if (power<=0) SETERRQ( PETSC_ERR_SUP_SYS, "Power must be 1 or up"); 848b2f1ab58SBarry Smith 849c328eeadSBarry Smith /* Copy input values*/ 850b2f1ab58SBarry Smith retval.nrows = ncols; 851b2f1ab58SBarry Smith retval.ncols = nrows; 852b2f1ab58SBarry Smith retval.nnz = 0; 853b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 854b2f1ab58SBarry Smith retval.block_data = PETSC_FALSE; 855b2f1ab58SBarry Smith 856c328eeadSBarry Smith /* Allocate sparseness pattern */ 85771d55d18SBarry Smith ierr = spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 858b2f1ab58SBarry Smith 859c328eeadSBarry Smith /* Allocate marker array */ 860b2f1ab58SBarry Smith ierr = PetscMalloc(nrows * sizeof(PetscInt), &iwork);CHKERRQ(ierr); 861b2f1ab58SBarry Smith 862c328eeadSBarry Smith /* Erase the pattern for this row */ 8634e1ff37aSBarry Smith ierr = PetscMemzero( (void *) iwork, retval.nrows*sizeof(PetscInt));CHKERRQ(ierr); 864b2f1ab58SBarry Smith 865c328eeadSBarry Smith /* Calculate marker values */ 866b2f1ab58SBarry Smith marker = 1; for (i=1; i<power; i++) {marker*=2;} 867b2f1ab58SBarry Smith 8686322f4bdSBarry Smith for (i=0; i<nrows; i++) { 869c328eeadSBarry Smith /* Calculate the pattern for each row */ 870b2f1ab58SBarry Smith 871b2f1ab58SBarry Smith nnz = in_matrix.row_nnz[i]; 872b2f1ab58SBarry Smith kend = i+in_matrix.icols[i][nnz-1]; 873b2f1ab58SBarry Smith if (maxmrk<=kend) {maxmrk=kend+1;} 87471d55d18SBarry Smith ierr = spbas_mark_row_power( iwork, i, &in_matrix, marker, i, maxmrk);CHKERRQ(ierr); 875b2f1ab58SBarry Smith 876c328eeadSBarry Smith /* Count the columns*/ 877b2f1ab58SBarry Smith nnz = 0; 878b2f1ab58SBarry Smith for (j=i; j<maxmrk; j++) {nnz+= (iwork[j]!=0);} 879b2f1ab58SBarry Smith 880c328eeadSBarry Smith /* Allocate the column indices */ 881b2f1ab58SBarry Smith retval.row_nnz[i] = nnz; 882b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&retval.icols[i]);CHKERRQ(ierr); 883b2f1ab58SBarry Smith 884c328eeadSBarry Smith /* Administrate the column indices */ 885b2f1ab58SBarry Smith inz = 0; 8866322f4bdSBarry Smith for (j=i; j<maxmrk; j++) { 8876322f4bdSBarry Smith if (iwork[j]) { 888b2f1ab58SBarry Smith retval.icols[i][inz] = j-i; 889b2f1ab58SBarry Smith inz++; 890b2f1ab58SBarry Smith iwork[j]=0; 891b2f1ab58SBarry Smith } 892b2f1ab58SBarry Smith } 893b2f1ab58SBarry Smith retval.nnz += nnz; 894b2f1ab58SBarry Smith }; 895b2f1ab58SBarry Smith ierr = PetscFree(iwork);CHKERRQ(ierr); 896b2f1ab58SBarry Smith *result = retval; 897b2f1ab58SBarry Smith PetscFunctionReturn(0); 898b2f1ab58SBarry Smith } 899b2f1ab58SBarry Smith 900b2f1ab58SBarry Smith 901b2f1ab58SBarry Smith 902b2f1ab58SBarry Smith /* 903b2f1ab58SBarry Smith spbas_keep_upper: 904b2f1ab58SBarry Smith remove the lower part of the matrix: keep the upper part 905b2f1ab58SBarry Smith */ 906b2f1ab58SBarry Smith #undef __FUNCT__ 907b2f1ab58SBarry Smith #define __FUNCT__ "spbas_keep_upper" 908b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper( spbas_matrix * inout_matrix) 909b2f1ab58SBarry Smith { 910b2f1ab58SBarry Smith PetscInt i, j; 911b2f1ab58SBarry Smith PetscInt jstart; 912b2f1ab58SBarry Smith 913b2f1ab58SBarry Smith PetscFunctionBegin; 9146322f4bdSBarry Smith if (inout_matrix->block_data) SETERRQ( PETSC_ERR_SUP_SYS, "Not yet for block data matrices\n"); 9156322f4bdSBarry Smith for (i=0; i<inout_matrix->nrows; i++) { 9166322f4bdSBarry Smith for (jstart=0; (jstart<inout_matrix->row_nnz[i]) && (inout_matrix->icols[i][jstart]<0); jstart++){} 9176322f4bdSBarry Smith if (jstart>0) { 9186322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9196322f4bdSBarry Smith inout_matrix->icols[i][j] = inout_matrix->icols[i][j+jstart]; 920b2f1ab58SBarry Smith } 921b2f1ab58SBarry Smith 9226322f4bdSBarry Smith if (inout_matrix->values) { 9236322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9246322f4bdSBarry Smith inout_matrix->values[i][j] = inout_matrix->values[i][j+jstart]; 925b2f1ab58SBarry Smith } 926b2f1ab58SBarry Smith } 927b2f1ab58SBarry Smith 928b2f1ab58SBarry Smith inout_matrix->row_nnz[i] -= jstart; 929b2f1ab58SBarry Smith 9306322f4bdSBarry Smith inout_matrix->icols[i] = (PetscInt*) realloc((void*) inout_matrix->icols[i], inout_matrix->row_nnz[i]*sizeof(PetscInt)); 931b2f1ab58SBarry Smith 9326322f4bdSBarry Smith if (inout_matrix->values) { 9336322f4bdSBarry Smith inout_matrix->values[i] = (PetscScalar*) realloc((void*) inout_matrix->values[i], inout_matrix->row_nnz[i]*sizeof(PetscScalar)); 934b2f1ab58SBarry Smith } 935b2f1ab58SBarry Smith inout_matrix->nnz -= jstart; 936b2f1ab58SBarry Smith } 937b2f1ab58SBarry Smith } 938b2f1ab58SBarry Smith PetscFunctionReturn(0); 939b2f1ab58SBarry Smith } 940b2f1ab58SBarry Smith 941