1b2f1ab58SBarry Smith #include "petsc.h" 2b2f1ab58SBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 3845006b9SBarry Smith #include "../src/mat/impls/aij/seq/bas/spbas.h" 4b2f1ab58SBarry Smith 5845006b9SBarry Smith /*MC 6845006b9SBarry Smith MAT_SOLVER_BAS - Provides ICC(k) with drop tolerance 7845006b9SBarry Smith 8845006b9SBarry Smith Works with MATAIJ matrices 9845006b9SBarry Smith 10845006b9SBarry Smith Options Database Keys: 11845006b9SBarry Smith + -pc_factor_levels <l> 12*b7c853c4SBarry Smith - -pc_factor_drop_tolerance 13845006b9SBarry Smith 14845006b9SBarry Smith Level: beginner 15845006b9SBarry Smith 16845006b9SBarry Smith Contributed by: Bas van 't Hof 17845006b9SBarry Smith 18*b7c853c4SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, PCFactorSetLevels(), PCFactorSetDropTolerance() 19845006b9SBarry Smith 20845006b9SBarry Smith M*/ 219ccc4a9bSBarry Smith 22b2f1ab58SBarry Smith /* 23b2f1ab58SBarry Smith spbas_memory_requirement: 24b2f1ab58SBarry Smith Calculate the number of bytes needed to store tha matrix 25b2f1ab58SBarry Smith */ 26b2f1ab58SBarry Smith #undef __FUNCT__ 27b2f1ab58SBarry Smith #define __FUNCT__ "spbas_memory_requirement" 28b2f1ab58SBarry Smith long int spbas_memory_requirement( spbas_matrix matrix) 29b2f1ab58SBarry Smith { 30c328eeadSBarry Smith long int memreq = 6 * sizeof(PetscInt) + /* nrows, ncols, nnz, n_alloc_icol, 31c328eeadSBarry Smith n_alloc_val, col_idx_type */ 32c328eeadSBarry Smith sizeof(PetscTruth) + /* block_data */ 33c328eeadSBarry Smith sizeof(PetscScalar**) + /* values */ 34c328eeadSBarry Smith sizeof(PetscScalar*) + /* alloc_val */ 35c328eeadSBarry Smith 2 * sizeof(int**) + /* icols, icols0 */ 36c328eeadSBarry Smith 2 * sizeof(PetscInt*) + /* row_nnz, alloc_icol */ 37c328eeadSBarry Smith matrix.nrows * sizeof(PetscInt) + /* row_nnz[*] */ 38c328eeadSBarry Smith matrix.nrows * sizeof(PetscInt*); /* icols[*] */ 39b2f1ab58SBarry Smith 40c328eeadSBarry Smith /* icol0[*] */ 414e1ff37aSBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) { memreq += matrix.nrows * sizeof(PetscInt); } 42b2f1ab58SBarry Smith 43c328eeadSBarry Smith /* icols[*][*] */ 444e1ff37aSBarry Smith if (matrix.block_data) { memreq += matrix.n_alloc_icol * sizeof(PetscInt); } 454e1ff37aSBarry Smith else { memreq += matrix.nnz * sizeof(PetscInt); } 46b2f1ab58SBarry Smith 474e1ff37aSBarry Smith if (matrix.values) { 48c328eeadSBarry Smith memreq += matrix.nrows * sizeof(PetscScalar*); /* values[*] */ 49c328eeadSBarry Smith /* values[*][*] */ 504e1ff37aSBarry Smith if (matrix.block_data) { memreq += matrix.n_alloc_val * sizeof(PetscScalar); } 514e1ff37aSBarry Smith else { memreq += matrix.nnz * sizeof(PetscScalar); } 52b2f1ab58SBarry Smith } 53b2f1ab58SBarry Smith return memreq; 54b2f1ab58SBarry Smith } 55b2f1ab58SBarry Smith 56b2f1ab58SBarry Smith /* 57b2f1ab58SBarry Smith spbas_allocate_pattern: 58b2f1ab58SBarry Smith allocate the pattern arrays row_nnz, icols and optionally values 59b2f1ab58SBarry Smith */ 60b2f1ab58SBarry Smith #undef __FUNCT__ 61b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_pattern" 62b2f1ab58SBarry Smith PetscErrorCode spbas_allocate_pattern( spbas_matrix * result, PetscTruth do_values) 63b2f1ab58SBarry Smith { 64b2f1ab58SBarry Smith PetscErrorCode ierr; 65b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 66b2f1ab58SBarry Smith PetscInt col_idx_type = result->col_idx_type; 67b2f1ab58SBarry Smith 68b2f1ab58SBarry Smith PetscFunctionBegin; 69c328eeadSBarry Smith /* Allocate sparseness pattern */ 70b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->row_nnz);CHKERRQ(ierr); 71b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt*),&result->icols);CHKERRQ(ierr); 72b2f1ab58SBarry Smith 73c328eeadSBarry Smith /* If offsets are given wrt an array, create array */ 744e1ff37aSBarry Smith if (col_idx_type == SPBAS_OFFSET_ARRAY) { 75b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->icol0);CHKERRQ(ierr); 764e1ff37aSBarry Smith } else { 77c328eeadSBarry Smith result->icol0 = PETSC_NULL; 78b2f1ab58SBarry Smith } 79b2f1ab58SBarry Smith 80c328eeadSBarry Smith /* If values are given, allocate values array */ 814e1ff37aSBarry Smith if (do_values) { 82c328eeadSBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscScalar*),&result->values);CHKERRQ(ierr); 834e1ff37aSBarry Smith } else { 84c328eeadSBarry Smith result->values = PETSC_NULL; 85b2f1ab58SBarry Smith } 86b2f1ab58SBarry Smith PetscFunctionReturn(0); 87b2f1ab58SBarry Smith } 88b2f1ab58SBarry Smith 89b2f1ab58SBarry Smith /* 90b2f1ab58SBarry Smith spbas_allocate_data: 91b2f1ab58SBarry Smith in case of block_data: 92b2f1ab58SBarry Smith Allocate the data arrays alloc_icol and optionally alloc_val, 93b2f1ab58SBarry Smith set appropriate pointers from icols and values; 94b2f1ab58SBarry Smith in case of !block_data: 95b2f1ab58SBarry Smith Allocate the arrays icols[i] and optionally values[i] 96b2f1ab58SBarry Smith */ 97b2f1ab58SBarry Smith #undef __FUNCT__ 98b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_data" 99b2f1ab58SBarry Smith PetscErrorCode spbas_allocate_data( spbas_matrix * result) 100b2f1ab58SBarry Smith { 101b2f1ab58SBarry Smith PetscInt i; 102b2f1ab58SBarry Smith PetscInt nnz = result->nnz; 103b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 104b2f1ab58SBarry Smith PetscInt r_nnz; 105b2f1ab58SBarry Smith PetscErrorCode ierr; 106c328eeadSBarry Smith PetscTruth do_values = (result->values != PETSC_NULL) ? PETSC_TRUE : PETSC_FALSE; 107b2f1ab58SBarry Smith PetscTruth block_data = result->block_data; 108b2f1ab58SBarry Smith 109b2f1ab58SBarry Smith PetscFunctionBegin; 1104e1ff37aSBarry Smith if (block_data) { 111c328eeadSBarry Smith /* Allocate the column number array and point to it */ 112b2f1ab58SBarry Smith result->n_alloc_icol = nnz; 113c328eeadSBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt), &result->alloc_icol);CHKERRQ(ierr); 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; 122c328eeadSBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscScalar), &result->alloc_val);CHKERRQ(ierr); 123b2f1ab58SBarry Smith result->values[0]= result->alloc_val; 1244e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 125b2f1ab58SBarry Smith result->values[i] = result->values[i-1] + result->row_nnz[i-1]; 126b2f1ab58SBarry Smith } 127b2f1ab58SBarry Smith } 1284e1ff37aSBarry Smith } else { 1294e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 130b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 131c328eeadSBarry Smith ierr = PetscMalloc(r_nnz*sizeof(PetscInt), &result->icols[i]);CHKERRQ(ierr); 132b2f1ab58SBarry Smith } 1334e1ff37aSBarry Smith if (do_values) { 1344e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 135b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 1364e1ff37aSBarry Smith ierr = PetscMalloc(r_nnz*sizeof(PetscScalar), &result->values[i]);CHKERRQ(ierr); 137b2f1ab58SBarry Smith } 138b2f1ab58SBarry Smith } 139b2f1ab58SBarry Smith } 140b2f1ab58SBarry Smith PetscFunctionReturn(0); 141b2f1ab58SBarry Smith } 142b2f1ab58SBarry Smith 143b2f1ab58SBarry Smith /* 144b2f1ab58SBarry Smith spbas_row_order_icol 145b2f1ab58SBarry Smith determine if row i1 should come 146b2f1ab58SBarry Smith + before row i2 in the sorted rows (return -1), 147b2f1ab58SBarry Smith + after (return 1) 148b2f1ab58SBarry Smith + is identical (return 0). 149b2f1ab58SBarry Smith */ 150b2f1ab58SBarry Smith #undef __FUNCT__ 151b2f1ab58SBarry Smith #define __FUNCT__ "spbas_row_order_icol" 152b2f1ab58SBarry Smith int spbas_row_order_icol(PetscInt i1, PetscInt i2, PetscInt *irow_in, PetscInt *icol_in,PetscInt col_idx_type) 153b2f1ab58SBarry Smith { 154b2f1ab58SBarry Smith PetscInt j; 155b2f1ab58SBarry Smith PetscInt nnz1 = irow_in[i1+1] - irow_in[i1]; 156b2f1ab58SBarry Smith PetscInt nnz2 = irow_in[i2+1] - irow_in[i2]; 157b2f1ab58SBarry Smith PetscInt * icol1 = &icol_in[irow_in[i1]]; 158b2f1ab58SBarry Smith PetscInt * icol2 = &icol_in[irow_in[i2]]; 159b2f1ab58SBarry Smith 160b2f1ab58SBarry Smith if (nnz1<nnz2) {return -1;} 161b2f1ab58SBarry Smith if (nnz1>nnz2) {return 1;} 162b2f1ab58SBarry Smith 1634e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 1644e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 165b2f1ab58SBarry Smith if (icol1[j]< icol2[j]) {return -1;} 166b2f1ab58SBarry Smith if (icol1[j]> icol2[j]) {return 1;} 167b2f1ab58SBarry Smith } 1684e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 1694e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 170b2f1ab58SBarry Smith if (icol1[j]-i1< icol2[j]-i2) {return -1;} 171b2f1ab58SBarry Smith if (icol1[j]-i1> icol2[j]-i2) {return 1;} 172b2f1ab58SBarry Smith } 1734e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 1744e1ff37aSBarry Smith for (j=1; j<nnz1; j++) { 175b2f1ab58SBarry Smith if (icol1[j]-icol1[0] < icol2[j]-icol2[0]) {return -1;} 176b2f1ab58SBarry Smith if (icol1[j]-icol1[0] > icol2[j]-icol2[0]) {return 1;} 177b2f1ab58SBarry Smith } 178b2f1ab58SBarry Smith } 179b2f1ab58SBarry Smith return 0; 180b2f1ab58SBarry Smith } 181b2f1ab58SBarry Smith 182b2f1ab58SBarry Smith /* 183b2f1ab58SBarry Smith spbas_mergesort_icols: 184b2f1ab58SBarry Smith return a sorting of the rows in which identical sparseness patterns are 185b2f1ab58SBarry Smith next to each other 186b2f1ab58SBarry Smith */ 187b2f1ab58SBarry Smith #undef __FUNCT__ 188b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort_icols" 189b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort_icols( PetscInt nrows, PetscInt * irow_in, PetscInt * icol_in,PetscInt col_idx_type, PetscInt *isort) 190b2f1ab58SBarry Smith { 191b2f1ab58SBarry Smith PetscErrorCode ierr; 192c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 193c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 194c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 195c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 196c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 197c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 198c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 199b2f1ab58SBarry Smith 200b2f1ab58SBarry Smith PetscFunctionBegin; 201b2f1ab58SBarry Smith 202b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&ialloc);CHKERRQ(ierr); 203b2f1ab58SBarry Smith 204b2f1ab58SBarry Smith ihlp1 = ialloc; 205b2f1ab58SBarry Smith ihlp2 = isort; 206b2f1ab58SBarry Smith 207c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 2084e1ff37aSBarry Smith for (istep=1; istep<nrows; istep*=2) { 209c328eeadSBarry Smith /* 210c328eeadSBarry Smith Combine sorted parts 211c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 212c328eeadSBarry Smith of ihlp2 and vhlp2 213c328eeadSBarry Smith 214c328eeadSBarry Smith into one sorted part 215c328eeadSBarry Smith istart:istart+2*istep-1 216c328eeadSBarry Smith of ihlp1 and vhlp1 217c328eeadSBarry Smith */ 2184e1ff37aSBarry Smith for (istart=0; istart<nrows; istart+=2*istep) { 219b2f1ab58SBarry Smith 220c328eeadSBarry Smith /* Set counters and bound array part endings */ 221b2f1ab58SBarry Smith i1=istart; i1end = i1+istep; if (i1end>nrows) {i1end=nrows;} 222b2f1ab58SBarry Smith 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) { 245b2f1ab58SBarry Smith 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; 271b2f1ab58SBarry Smith PetscTruth *used; 272b2f1ab58SBarry Smith PetscInt ptr; 273b2f1ab58SBarry Smith PetscInt i,j; 274b2f1ab58SBarry Smith const PetscTruth 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; 283b2f1ab58SBarry Smith ierr = spbas_allocate_pattern( B, no_values);CHKERRQ(ierr); 284b2f1ab58SBarry Smith 285c328eeadSBarry Smith /* When using an offset array, set it */ 2864e1ff37aSBarry Smith if (col_idx_type==SPBAS_OFFSET_ARRAY) { 287b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {B->icol0[i] = icol_in[irow_in[i]];} 288b2f1ab58SBarry Smith } 289b2f1ab58SBarry Smith 290c328eeadSBarry Smith /* Allocate the ordering for the rows */ 291b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&isort);CHKERRQ(ierr); 292b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&ipoint);CHKERRQ(ierr); 293b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscTruth),&used);CHKERRQ(ierr); 294b2f1ab58SBarry Smith 295c328eeadSBarry Smith /* Initialize the sorting */ 2964e1ff37aSBarry Smith ierr = PetscMemzero((void*) used, nrows*sizeof(PetscTruth));CHKERRQ(ierr); 2974e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 298b2f1ab58SBarry Smith B->row_nnz[i] = irow_in[i+1]-irow_in[i]; 299b2f1ab58SBarry Smith isort[i] = i; 300b2f1ab58SBarry Smith ipoint[i]= i; 301b2f1ab58SBarry Smith } 302b2f1ab58SBarry Smith 303c328eeadSBarry Smith /* Sort the rows so that identical columns will be next to each other */ 304b2f1ab58SBarry Smith ierr = spbas_mergesort_icols( nrows, irow_in, icol_in, col_idx_type, isort);CHKERRQ(ierr); 30571d55d18SBarry Smith ierr = PetscInfo(PETSC_NULL,"Rows have been sorted for patterns\n");CHKERRQ(ierr); 306b2f1ab58SBarry Smith 307c328eeadSBarry Smith /* Replace identical rows with the first one in the list */ 3084e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 3094e1ff37aSBarry Smith if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) { 310b2f1ab58SBarry Smith ipoint[isort[i]] = ipoint[isort[i-1]]; 311b2f1ab58SBarry Smith } 312b2f1ab58SBarry Smith } 313b2f1ab58SBarry Smith 314c328eeadSBarry Smith /* Collect the rows which are used*/ 315b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {used[ipoint[i]] = PETSC_TRUE;} 316b2f1ab58SBarry Smith 317c328eeadSBarry Smith /* Calculate needed memory */ 318b2f1ab58SBarry Smith B->n_alloc_icol = 0; 3194e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 320b2f1ab58SBarry Smith if (used[i]) {B->n_alloc_icol += B->row_nnz[i];} 321b2f1ab58SBarry Smith } 32271d55d18SBarry Smith ierr = PetscMalloc(B->n_alloc_icol*sizeof(PetscInt),&B->alloc_icol);CHKERRQ(ierr); 323b2f1ab58SBarry Smith 324c328eeadSBarry Smith /* Fill in the diagonal offsets for the rows which store their own data */ 325b2f1ab58SBarry Smith ptr = 0; 3264e1ff37aSBarry Smith for (i=0; i<B->nrows; i++) { 3274e1ff37aSBarry Smith if (used[i]) { 328b2f1ab58SBarry Smith B->icols[i] = &B->alloc_icol[ptr]; 329b2f1ab58SBarry Smith icols = &icol_in[irow_in[i]]; 330b2f1ab58SBarry Smith row_nnz = B->row_nnz[i]; 3314e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 3324e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 333b2f1ab58SBarry Smith B->icols[i][j] = icols[j]; 334b2f1ab58SBarry Smith } 3354e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 3364e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 337b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-i; 338b2f1ab58SBarry Smith } 3394e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 3404e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 341b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-icols[0]; 342b2f1ab58SBarry Smith } 343b2f1ab58SBarry Smith } 344b2f1ab58SBarry Smith ptr += B->row_nnz[i]; 345b2f1ab58SBarry Smith } 346b2f1ab58SBarry Smith } 347b2f1ab58SBarry Smith 348c328eeadSBarry Smith /* Point to the right places for all data */ 3494e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 350b2f1ab58SBarry Smith B->icols[i] = B->icols[ipoint[i]]; 351b2f1ab58SBarry Smith } 352c328eeadSBarry Smith ierr = PetscInfo(PETSC_NULL,"Row patterns have been compressed\n");CHKERRQ(ierr); 353c328eeadSBarry Smith ierr = PetscInfo1(PETSC_NULL," (%G nonzeros per row)\n", (PetscReal) nnz / (PetscReal) nrows);CHKERRQ(ierr); 354b2f1ab58SBarry Smith 355b2f1ab58SBarry Smith ierr=PetscFree(isort);CHKERRQ(ierr); 356b2f1ab58SBarry Smith ierr=PetscFree(used);CHKERRQ(ierr); 357b2f1ab58SBarry Smith ierr=PetscFree(ipoint);CHKERRQ(ierr); 358b2f1ab58SBarry Smith 359b2f1ab58SBarry Smith mem_compressed = spbas_memory_requirement( *B ); 3604e1ff37aSBarry Smith *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig; 361b2f1ab58SBarry Smith PetscFunctionReturn(0); 362b2f1ab58SBarry Smith } 363b2f1ab58SBarry Smith 364b2f1ab58SBarry Smith /* 365b2f1ab58SBarry Smith spbas_incomplete_cholesky 366b2f1ab58SBarry Smith Incomplete Cholesky decomposition 367b2f1ab58SBarry Smith */ 368845006b9SBarry Smith #include "../src/mat/impls/aij/seq/bas/spbas_cholesky.h" 369b2f1ab58SBarry Smith 370b2f1ab58SBarry Smith /* 371b2f1ab58SBarry Smith spbas_delete : de-allocate the arrays owned by this matrix 372b2f1ab58SBarry Smith */ 373b2f1ab58SBarry Smith #undef __FUNCT__ 374b2f1ab58SBarry Smith #define __FUNCT__ "spbas_delete" 375b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix matrix) 376b2f1ab58SBarry Smith { 377b2f1ab58SBarry Smith PetscInt i; 378b2f1ab58SBarry Smith PetscErrorCode ierr; 3799ccc4a9bSBarry Smith 380b2f1ab58SBarry Smith PetscFunctionBegin; 3819ccc4a9bSBarry Smith if (matrix.block_data) { 382b2f1ab58SBarry Smith ierr=PetscFree(matrix.alloc_icol);CHKERRQ(ierr) 383b2f1ab58SBarry Smith if (matrix.values){ierr=PetscFree(matrix.alloc_val);CHKERRQ(ierr);} 3849ccc4a9bSBarry Smith } else { 3859ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.icols[i]);CHKERRQ(ierr);} 386b2f1ab58SBarry Smith ierr = PetscFree(matrix.icols);CHKERRQ(ierr); 3879ccc4a9bSBarry Smith if (matrix.values) { 3889ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.values[i]);CHKERRQ(ierr);} 389b2f1ab58SBarry Smith } 390b2f1ab58SBarry Smith } 391b2f1ab58SBarry Smith 392b2f1ab58SBarry Smith ierr=PetscFree(matrix.row_nnz);CHKERRQ(ierr); 393b2f1ab58SBarry Smith ierr=PetscFree(matrix.icols);CHKERRQ(ierr); 3949ccc4a9bSBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) {ierr=PetscFree(matrix.icol0);CHKERRQ(ierr);} 3959ccc4a9bSBarry Smith if (matrix.values) {ierr=PetscFree(matrix.values);CHKERRQ(ierr);} 396b2f1ab58SBarry Smith PetscFunctionReturn(0); 397b2f1ab58SBarry Smith } 398b2f1ab58SBarry Smith 399b2f1ab58SBarry Smith /* 400b2f1ab58SBarry Smith spbas_matrix_to_crs: 401b2f1ab58SBarry Smith Convert an spbas_matrix to compessed row storage 402b2f1ab58SBarry Smith */ 403b2f1ab58SBarry Smith #undef __FUNCT__ 404b2f1ab58SBarry Smith #define __FUNCT__ "spbas_matrix_to_crs" 405b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out) 406b2f1ab58SBarry Smith { 407b2f1ab58SBarry Smith PetscInt nrows = matrix_A.nrows; 408b2f1ab58SBarry Smith PetscInt nnz = matrix_A.nnz; 409b2f1ab58SBarry Smith PetscInt i,j,r_nnz,i0; 410b2f1ab58SBarry Smith PetscInt *irow; 411b2f1ab58SBarry Smith PetscInt *icol; 412b2f1ab58SBarry Smith PetscInt *icol_A; 413b2f1ab58SBarry Smith MatScalar *val; 414b2f1ab58SBarry Smith PetscScalar *val_A; 415b2f1ab58SBarry Smith PetscInt col_idx_type = matrix_A.col_idx_type; 416d36aa34eSBarry Smith PetscTruth do_values = matrix_A.values ? PETSC_TRUE : PETSC_FALSE; 417b2f1ab58SBarry Smith PetscErrorCode ierr; 418b2f1ab58SBarry Smith 419b2f1ab58SBarry Smith PetscFunctionBegin; 420b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt) * (nrows+1), &irow);CHKERRQ(ierr); 421b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt) * nnz, &icol);CHKERRQ(ierr); 4229ccc4a9bSBarry Smith *icol_out = icol; 4239ccc4a9bSBarry Smith *irow_out=irow; 4249ccc4a9bSBarry Smith if (do_values) { 425b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(MatScalar) * nnz, &val);CHKERRQ(ierr); 426b2f1ab58SBarry Smith *val_out = val; *icol_out = icol; *irow_out=irow; 427b2f1ab58SBarry Smith } 428b2f1ab58SBarry Smith 429b2f1ab58SBarry Smith irow[0]=0; 4309ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 431b2f1ab58SBarry Smith r_nnz = matrix_A.row_nnz[i]; 432b2f1ab58SBarry Smith i0 = irow[i]; 433b2f1ab58SBarry Smith irow[i+1] = i0 + r_nnz; 434b2f1ab58SBarry Smith icol_A = matrix_A.icols[i]; 435b2f1ab58SBarry Smith 4369ccc4a9bSBarry Smith if (do_values) { 437b2f1ab58SBarry Smith val_A = matrix_A.values[i]; 4389ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 439b2f1ab58SBarry Smith icol[i0+j] = icol_A[j]; 440b2f1ab58SBarry Smith val[i0+j] = val_A[j]; 441b2f1ab58SBarry Smith } 4429ccc4a9bSBarry Smith } else { 443b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] = icol_A[j]; } 444b2f1ab58SBarry Smith } 445b2f1ab58SBarry Smith 4469ccc4a9bSBarry Smith if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 447b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] += i; } 448b2f1ab58SBarry Smith } 4499ccc4a9bSBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 450b2f1ab58SBarry Smith i0 = matrix_A.icol0[i]; 451b2f1ab58SBarry Smith 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; 473b2f1ab58SBarry Smith PetscInt icol0; 474b2f1ab58SBarry Smith PetscScalar * val; 475b2f1ab58SBarry Smith PetscErrorCode ierr; 4764e1ff37aSBarry Smith 477b2f1ab58SBarry Smith PetscFunctionBegin; 478b2f1ab58SBarry Smith 479c328eeadSBarry Smith /* Copy input values */ 480b2f1ab58SBarry Smith result->nrows = nrows; 481b2f1ab58SBarry Smith result->ncols = ncols; 482b2f1ab58SBarry Smith result->nnz = nnz; 483b2f1ab58SBarry Smith result->col_idx_type = SPBAS_COLUMN_NUMBERS; 484b2f1ab58SBarry Smith result->block_data = PETSC_TRUE; 485b2f1ab58SBarry Smith 486c328eeadSBarry Smith /* Allocate sparseness pattern */ 48771d55d18SBarry Smith ierr = spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 488b2f1ab58SBarry Smith 489c328eeadSBarry Smith /* Count the number of nonzeros in each row */ 490b2f1ab58SBarry Smith for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; } 491b2f1ab58SBarry Smith 4929ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 493b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 494b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 4959ccc4a9bSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 496b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[irow[j]]++; } 4979ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 498b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[i+irow[j]]++; } 4999ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 500b2f1ab58SBarry Smith icol0=in_matrix.icol0[i]; 501b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[icol0+irow[j]]++; } 502b2f1ab58SBarry Smith } 503b2f1ab58SBarry Smith } 504b2f1ab58SBarry Smith 505c328eeadSBarry Smith /* Set the pointers to the data */ 506b2f1ab58SBarry Smith ierr = spbas_allocate_data(result);CHKERRQ(ierr); 507b2f1ab58SBarry Smith 508c328eeadSBarry Smith /* Reset the number of nonzeros in each row */ 509b2f1ab58SBarry Smith for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; } 510b2f1ab58SBarry Smith 511c328eeadSBarry Smith /* Fill the data arrays */ 5129ccc4a9bSBarry Smith if (in_matrix.values) { 5139ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 514b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 515b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 516b2f1ab58SBarry Smith val = in_matrix.values[i]; 517b2f1ab58SBarry Smith 518b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) {icol0=0;} 519b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;} 520b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) 521b2f1ab58SBarry Smith {icol0=in_matrix.icol0[i];} 5229ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 523b2f1ab58SBarry Smith k = icol0 + irow[j]; 524b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 525b2f1ab58SBarry Smith result->values[k][result->row_nnz[k]] = val[j]; 526b2f1ab58SBarry Smith result->row_nnz[k]++; 527b2f1ab58SBarry Smith } 528b2f1ab58SBarry Smith } 5299ccc4a9bSBarry Smith } else { 5309ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 531b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 532b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 533b2f1ab58SBarry Smith 534b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) {icol0=0;} 535b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;} 5369ccc4a9bSBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) {icol0=in_matrix.icol0[i];} 537b2f1ab58SBarry Smith 5389ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 539b2f1ab58SBarry Smith k = icol0 + irow[j]; 540b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 541b2f1ab58SBarry Smith result->row_nnz[k]++; 542b2f1ab58SBarry Smith } 543b2f1ab58SBarry Smith } 544b2f1ab58SBarry Smith } 545b2f1ab58SBarry Smith PetscFunctionReturn(0); 546b2f1ab58SBarry Smith } 547b2f1ab58SBarry Smith 548b2f1ab58SBarry Smith /* 549b2f1ab58SBarry Smith spbas_mergesort 550b2f1ab58SBarry Smith 551b2f1ab58SBarry Smith mergesort for an array of intergers and an array of associated 552b2f1ab58SBarry Smith reals 553b2f1ab58SBarry Smith 554b2f1ab58SBarry Smith on output, icol[0..nnz-1] is increasing; 555b2f1ab58SBarry Smith val[0..nnz-1] has undergone the same permutation as icol 556b2f1ab58SBarry Smith 557c328eeadSBarry Smith NB: val may be PETSC_NULL: in that case, only the integers are sorted 558b2f1ab58SBarry Smith 559b2f1ab58SBarry Smith */ 560b2f1ab58SBarry Smith #undef __FUNCT__ 561b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort" 562b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val) 563b2f1ab58SBarry Smith { 564c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 565c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 566c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 567c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 568c328eeadSBarry Smith PetscScalar *valloc=PETSC_NULL; 569c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 570b2f1ab58SBarry Smith PetscScalar *vswap; 571c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 572c328eeadSBarry Smith PetscScalar *vhlp1=PETSC_NULL; /* (arrays under construction) */ 573c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 574c328eeadSBarry Smith PetscScalar *vhlp2=PETSC_NULL; 575b2f1ab58SBarry Smith PetscErrorCode ierr; 576b2f1ab58SBarry Smith 577b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&ialloc);CHKERRQ(ierr); 578b2f1ab58SBarry Smith ihlp1 = ialloc; 579b2f1ab58SBarry Smith ihlp2 = icol; 580b2f1ab58SBarry Smith 5819ccc4a9bSBarry Smith if (val) { 582b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscScalar),&valloc);CHKERRQ(ierr); 583b2f1ab58SBarry Smith vhlp1 = valloc; 584b2f1ab58SBarry Smith vhlp2 = val; 585b2f1ab58SBarry Smith } 586b2f1ab58SBarry Smith 587b2f1ab58SBarry Smith 588c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 5899ccc4a9bSBarry Smith for (istep=1; istep<nnz; istep*=2) { 590c328eeadSBarry Smith /* 591c328eeadSBarry Smith Combine sorted parts 592c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 593c328eeadSBarry Smith of ihlp2 and vhlp2 594c328eeadSBarry Smith 595c328eeadSBarry Smith into one sorted part 596c328eeadSBarry Smith istart:istart+2*istep-1 597c328eeadSBarry Smith of ihlp1 and vhlp1 598c328eeadSBarry Smith */ 5999ccc4a9bSBarry Smith for (istart=0; istart<nnz; istart+=2*istep) { 600b2f1ab58SBarry Smith 601c328eeadSBarry Smith /* Set counters and bound array part endings */ 602b2f1ab58SBarry Smith i1=istart; i1end = i1+istep; if (i1end>nnz) {i1end=nnz;} 603b2f1ab58SBarry Smith i2=istart+istep; i2end = i2+istep; if (i2end>nnz) {i2end=nnz;} 604b2f1ab58SBarry Smith 605c328eeadSBarry Smith /* Merge the two array parts */ 6069ccc4a9bSBarry Smith if (val) { 6079ccc4a9bSBarry Smith for (i=istart; i<i2end; i++) { 6089ccc4a9bSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 609b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 610b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 611b2f1ab58SBarry Smith i1++; 6129ccc4a9bSBarry Smith } else if (i2<i2end ){ 613b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 614b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i2]; 615b2f1ab58SBarry Smith i2++; 6169ccc4a9bSBarry Smith } else { 617b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 618b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 619b2f1ab58SBarry Smith i1++; 620b2f1ab58SBarry Smith } 621b2f1ab58SBarry Smith } 6229ccc4a9bSBarry Smith } else { 6236322f4bdSBarry Smith for (i=istart; i<i2end; i++) { 6246322f4bdSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 625b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 626b2f1ab58SBarry Smith i1++; 6276322f4bdSBarry Smith } else if (i2<i2end ) { 628b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 629b2f1ab58SBarry Smith i2++; 6306322f4bdSBarry Smith } else { 631b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 632b2f1ab58SBarry Smith i1++; 633b2f1ab58SBarry Smith } 634b2f1ab58SBarry Smith } 635b2f1ab58SBarry Smith } 636b2f1ab58SBarry Smith } 637b2f1ab58SBarry Smith 638c328eeadSBarry Smith /* Swap the two array sets */ 639b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 640b2f1ab58SBarry Smith vswap = vhlp2; vhlp2 = vhlp1; vhlp1 = vswap; 641b2f1ab58SBarry Smith } 642b2f1ab58SBarry Smith 643c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 6446322f4bdSBarry Smith if (ihlp2 != icol) { 645b2f1ab58SBarry Smith for (i=0; i<nnz; i++) { icol[i] = ihlp2[i]; } 6466322f4bdSBarry Smith if (val) { 647b2f1ab58SBarry Smith for (i=0; i<nnz; i++) { val[i] = vhlp2[i]; } 648b2f1ab58SBarry Smith } 649b2f1ab58SBarry Smith } 650b2f1ab58SBarry Smith 651b2f1ab58SBarry Smith ierr = PetscFree(ialloc);CHKERRQ(ierr); 652b2f1ab58SBarry Smith if(val){ierr = PetscFree(valloc);CHKERRQ(ierr);} 653b2f1ab58SBarry Smith PetscFunctionReturn(0); 654b2f1ab58SBarry Smith } 655b2f1ab58SBarry Smith 656b2f1ab58SBarry Smith /* 657b2f1ab58SBarry Smith spbas_apply_reordering_rows: 658b2f1ab58SBarry Smith apply the given reordering to the rows: matrix_A = matrix_A(perm,:); 659b2f1ab58SBarry Smith */ 660b2f1ab58SBarry Smith #undef __FUNCT__ 661b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_rows" 662b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_rows(spbas_matrix *matrix_A, const PetscInt *permutation) 663b2f1ab58SBarry Smith { 664b2f1ab58SBarry Smith PetscInt i,j,ip; 665b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 666b2f1ab58SBarry Smith PetscInt * row_nnz; 667b2f1ab58SBarry Smith PetscInt **icols; 668d36aa34eSBarry Smith PetscTruth do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 669c328eeadSBarry Smith PetscScalar **vals=PETSC_NULL; 670b2f1ab58SBarry Smith PetscErrorCode ierr; 671b2f1ab58SBarry Smith 672b2f1ab58SBarry Smith PetscFunctionBegin; 6736322f4bdSBarry Smith if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ( PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern\n"); 674b2f1ab58SBarry Smith 6756322f4bdSBarry Smith if (do_values) { 676b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscScalar*)*nrows, &vals);CHKERRQ(ierr); 677b2f1ab58SBarry Smith } 678b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt)*nrows, &row_nnz);CHKERRQ(ierr); 679b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt*)*nrows, &icols);CHKERRQ(ierr); 680b2f1ab58SBarry Smith 6816322f4bdSBarry Smith for (i=0; i<nrows;i++) { 682b2f1ab58SBarry Smith ip = permutation[i]; 683b2f1ab58SBarry Smith if (do_values) {vals[i] = matrix_A->values[ip];} 684b2f1ab58SBarry Smith icols[i] = matrix_A->icols[ip]; 685b2f1ab58SBarry Smith row_nnz[i] = matrix_A->row_nnz[ip]; 686b2f1ab58SBarry Smith for (j=0; j<row_nnz[i]; j++) { icols[i][j] += ip-i; } 687b2f1ab58SBarry Smith } 688b2f1ab58SBarry Smith 689b2f1ab58SBarry Smith if (do_values){ ierr = PetscFree(matrix_A->values);CHKERRQ(ierr);} 690b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->icols);CHKERRQ(ierr); 691b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->row_nnz);CHKERRQ(ierr); 692b2f1ab58SBarry Smith 693b2f1ab58SBarry Smith if (do_values) { matrix_A->values = vals; } 694b2f1ab58SBarry Smith matrix_A->icols = icols; 695b2f1ab58SBarry Smith matrix_A->row_nnz = row_nnz; 696b2f1ab58SBarry Smith 697b2f1ab58SBarry Smith PetscFunctionReturn(0); 698b2f1ab58SBarry Smith } 699b2f1ab58SBarry Smith 700b2f1ab58SBarry Smith 701b2f1ab58SBarry Smith /* 702b2f1ab58SBarry Smith spbas_apply_reordering_cols: 703b2f1ab58SBarry Smith apply the given reordering to the columns: matrix_A(:,perm) = matrix_A; 704b2f1ab58SBarry Smith */ 705b2f1ab58SBarry Smith #undef __FUNCT__ 706b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_cols" 707b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_cols( spbas_matrix *matrix_A,const PetscInt *permutation) 708b2f1ab58SBarry Smith { 709b2f1ab58SBarry Smith PetscInt i,j; 710b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 711b2f1ab58SBarry Smith PetscInt row_nnz; 712b2f1ab58SBarry Smith PetscInt *icols; 713d36aa34eSBarry Smith PetscTruth do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 714c328eeadSBarry Smith PetscScalar *vals=PETSC_NULL; 715b2f1ab58SBarry Smith PetscErrorCode ierr; 716b2f1ab58SBarry Smith 717b2f1ab58SBarry Smith PetscFunctionBegin; 7186322f4bdSBarry Smith if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ( PETSC_ERR_SUP_SYS, "must have diagonal offsets in pattern\n"); 719b2f1ab58SBarry Smith 7206322f4bdSBarry Smith for (i=0; i<nrows;i++) { 721b2f1ab58SBarry Smith icols = matrix_A->icols[i]; 722b2f1ab58SBarry Smith row_nnz = matrix_A->row_nnz[i]; 723b2f1ab58SBarry Smith if (do_values) { vals = matrix_A->values[i]; } 724b2f1ab58SBarry Smith 7256322f4bdSBarry Smith for (j=0; j<row_nnz; j++) { 726b2f1ab58SBarry Smith icols[j] = permutation[i+icols[j]]-i; 727b2f1ab58SBarry Smith } 728b2f1ab58SBarry Smith ierr = spbas_mergesort(row_nnz, icols, vals);CHKERRQ(ierr); 729b2f1ab58SBarry Smith } 730b2f1ab58SBarry Smith PetscFunctionReturn(0); 731b2f1ab58SBarry Smith } 732b2f1ab58SBarry Smith 733b2f1ab58SBarry Smith /* 734b2f1ab58SBarry Smith spbas_apply_reordering: 735b2f1ab58SBarry Smith apply the given reordering: matrix_A(perm,perm) = matrix_A; 736b2f1ab58SBarry Smith */ 737b2f1ab58SBarry Smith #undef __FUNCT__ 738b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering" 739b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering( spbas_matrix *matrix_A, const PetscInt *permutation, const PetscInt * inv_perm) 740b2f1ab58SBarry Smith { 741b2f1ab58SBarry Smith PetscErrorCode ierr; 742b2f1ab58SBarry Smith PetscFunctionBegin; 743b2f1ab58SBarry Smith ierr = spbas_apply_reordering_rows( matrix_A, inv_perm);CHKERRQ(ierr); 744b2f1ab58SBarry Smith ierr = spbas_apply_reordering_cols( matrix_A, permutation);CHKERRQ(ierr); 745b2f1ab58SBarry Smith PetscFunctionReturn(0); 746b2f1ab58SBarry Smith } 747b2f1ab58SBarry Smith 748b2f1ab58SBarry Smith #undef __FUNCT__ 749b2f1ab58SBarry Smith #define __FUNCT__ "spbas_pattern_only" 750b2f1ab58SBarry Smith PetscErrorCode spbas_pattern_only(PetscInt nrows, PetscInt ncols, PetscInt *ai, PetscInt *aj, spbas_matrix * result) 751b2f1ab58SBarry Smith { 752b2f1ab58SBarry Smith spbas_matrix retval; 753b2f1ab58SBarry Smith PetscInt i, j, i0, r_nnz; 754b2f1ab58SBarry Smith PetscErrorCode ierr; 755b2f1ab58SBarry Smith 756b2f1ab58SBarry Smith PetscFunctionBegin; 757b2f1ab58SBarry Smith 758c328eeadSBarry Smith /* Copy input values */ 759b2f1ab58SBarry Smith retval.nrows = nrows; 760b2f1ab58SBarry Smith retval.ncols = ncols; 761b2f1ab58SBarry Smith retval.nnz = ai[nrows]; 762b2f1ab58SBarry Smith 763b2f1ab58SBarry Smith retval.block_data = PETSC_TRUE; 764b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 765b2f1ab58SBarry Smith 766c328eeadSBarry Smith /* Allocate output matrix */ 767b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(&retval, PETSC_FALSE);CHKERRQ(ierr); 768b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {retval.row_nnz[i] = ai[i+1]-ai[i];} 769b2f1ab58SBarry Smith ierr = spbas_allocate_data(&retval);CHKERRQ(ierr); 770c328eeadSBarry Smith /* Copy the structure */ 7716322f4bdSBarry Smith for (i = 0; i<retval.nrows; i++) { 772b2f1ab58SBarry Smith i0 = ai[i]; 773b2f1ab58SBarry Smith r_nnz = ai[i+1]-i0; 774b2f1ab58SBarry Smith 7756322f4bdSBarry Smith for (j=0; j<r_nnz; j++) { 776b2f1ab58SBarry Smith retval.icols[i][j] = aj[i0+j]-i; 777b2f1ab58SBarry Smith } 778b2f1ab58SBarry Smith } 779b2f1ab58SBarry Smith *result = retval; 780b2f1ab58SBarry Smith PetscFunctionReturn(0); 781b2f1ab58SBarry Smith } 782b2f1ab58SBarry Smith 783b2f1ab58SBarry Smith 784b2f1ab58SBarry Smith /* 785b2f1ab58SBarry Smith spbas_mark_row_power: 786b2f1ab58SBarry Smith Mark the columns in row 'row' which are nonzero in 787b2f1ab58SBarry Smith matrix^2log(marker). 788b2f1ab58SBarry Smith */ 789b2f1ab58SBarry Smith #undef __FUNCT__ 790b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mark_row_power" 791b2f1ab58SBarry Smith PetscErrorCode spbas_mark_row_power( 792c328eeadSBarry Smith PetscInt *iwork, /* marker-vector */ 793c328eeadSBarry Smith PetscInt row, /* row for which the columns are marked */ 794c328eeadSBarry Smith spbas_matrix * in_matrix, /* matrix for which the power is being calculated */ 795c328eeadSBarry Smith PetscInt marker, /* marker-value: 2^power */ 796c328eeadSBarry Smith PetscInt minmrk, /* lower bound for marked points */ 797c328eeadSBarry Smith PetscInt maxmrk) /* upper bound for marked points */ 798b2f1ab58SBarry Smith { 799b2f1ab58SBarry Smith PetscErrorCode ierr; 800b2f1ab58SBarry Smith PetscInt i,j, nnz; 801b2f1ab58SBarry Smith 802b2f1ab58SBarry Smith PetscFunctionBegin; 803b2f1ab58SBarry Smith nnz = in_matrix->row_nnz[row]; 804b2f1ab58SBarry Smith 805c328eeadSBarry Smith /* For higher powers, call this function recursively */ 8066322f4bdSBarry Smith if (marker>1) { 8076322f4bdSBarry Smith for (i=0; i<nnz;i++) { 808b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 8096322f4bdSBarry Smith if (minmrk<=j && j<maxmrk && iwork[j] < marker ) { 81071d55d18SBarry Smith ierr = spbas_mark_row_power( iwork, row + in_matrix->icols[row][i],in_matrix, marker/2,minmrk,maxmrk);CHKERRQ(ierr); 811b2f1ab58SBarry Smith iwork[j] |= marker; 812b2f1ab58SBarry Smith } 813b2f1ab58SBarry Smith } 8146322f4bdSBarry Smith } else { 815c328eeadSBarry Smith /* Mark the columns reached */ 8166322f4bdSBarry Smith for (i=0; i<nnz;i++) { 817b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 818b2f1ab58SBarry Smith if (minmrk<=j && j<maxmrk) { iwork[j] |= 1; } 819b2f1ab58SBarry Smith } 820b2f1ab58SBarry Smith } 821b2f1ab58SBarry Smith PetscFunctionReturn(0); 822b2f1ab58SBarry Smith } 823b2f1ab58SBarry Smith 824b2f1ab58SBarry Smith 825b2f1ab58SBarry Smith /* 826b2f1ab58SBarry Smith spbas_power 827b2f1ab58SBarry Smith Calculate sparseness patterns for incomplete Cholesky decompositions 828b2f1ab58SBarry Smith of a given order: (almost) all nonzeros of the matrix^(order+1) which 829b2f1ab58SBarry Smith are inside the band width are found and stored in the output sparseness 830b2f1ab58SBarry Smith pattern. 831b2f1ab58SBarry Smith */ 832b2f1ab58SBarry Smith #undef __FUNCT__ 833b2f1ab58SBarry Smith #define __FUNCT__ "spbas_power" 834b2f1ab58SBarry Smith PetscErrorCode spbas_power (spbas_matrix in_matrix,PetscInt power, spbas_matrix * result) 835b2f1ab58SBarry Smith { 836b2f1ab58SBarry Smith spbas_matrix retval; 837b2f1ab58SBarry Smith PetscInt nrows = in_matrix.nrows; 838b2f1ab58SBarry Smith PetscInt ncols = in_matrix.ncols; 839b2f1ab58SBarry Smith PetscInt i, j, kend; 840b2f1ab58SBarry Smith PetscInt nnz, inz; 841b2f1ab58SBarry Smith PetscInt *iwork; 842b2f1ab58SBarry Smith PetscInt marker; 843b2f1ab58SBarry Smith PetscInt maxmrk=0; 844b2f1ab58SBarry Smith PetscErrorCode ierr; 845b2f1ab58SBarry Smith 846b2f1ab58SBarry Smith PetscFunctionBegin; 8476322f4bdSBarry Smith if (in_matrix.col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ( PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern\n"); 8486322f4bdSBarry Smith if (ncols != nrows) SETERRQ( PETSC_ERR_ARG_INCOMP, "Dimension error\n"); 8496322f4bdSBarry Smith if (in_matrix.values) SETERRQ( PETSC_ERR_ARG_INCOMP, "Input array must be sparseness pattern (no values)"); 8506322f4bdSBarry Smith if (power<=0) SETERRQ( PETSC_ERR_SUP_SYS, "Power must be 1 or up"); 851b2f1ab58SBarry Smith 852c328eeadSBarry Smith /* Copy input values*/ 853b2f1ab58SBarry Smith retval.nrows = ncols; 854b2f1ab58SBarry Smith retval.ncols = nrows; 855b2f1ab58SBarry Smith retval.nnz = 0; 856b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 857b2f1ab58SBarry Smith retval.block_data = PETSC_FALSE; 858b2f1ab58SBarry Smith 859c328eeadSBarry Smith /* Allocate sparseness pattern */ 86071d55d18SBarry Smith ierr = spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 861b2f1ab58SBarry Smith 862c328eeadSBarry Smith /* Allocate marker array */ 863b2f1ab58SBarry Smith ierr = PetscMalloc(nrows * sizeof(PetscInt), &iwork);CHKERRQ(ierr); 864b2f1ab58SBarry Smith 865c328eeadSBarry Smith /* Erase the pattern for this row */ 8664e1ff37aSBarry Smith ierr = PetscMemzero( (void *) iwork, retval.nrows*sizeof(PetscInt));CHKERRQ(ierr); 867b2f1ab58SBarry Smith 868c328eeadSBarry Smith /* Calculate marker values */ 869b2f1ab58SBarry Smith marker = 1; for (i=1; i<power; i++) {marker*=2;} 870b2f1ab58SBarry Smith 8716322f4bdSBarry Smith for (i=0; i<nrows; i++) { 872c328eeadSBarry Smith /* Calculate the pattern for each row */ 873b2f1ab58SBarry Smith 874b2f1ab58SBarry Smith nnz = in_matrix.row_nnz[i]; 875b2f1ab58SBarry Smith kend = i+in_matrix.icols[i][nnz-1]; 876b2f1ab58SBarry Smith if (maxmrk<=kend) {maxmrk=kend+1;} 87771d55d18SBarry Smith ierr = spbas_mark_row_power( iwork, i, &in_matrix, marker, i, maxmrk);CHKERRQ(ierr); 878b2f1ab58SBarry Smith 879c328eeadSBarry Smith /* Count the columns*/ 880b2f1ab58SBarry Smith nnz = 0; 881b2f1ab58SBarry Smith for (j=i; j<maxmrk; j++) {nnz+= (iwork[j]!=0);} 882b2f1ab58SBarry Smith 883c328eeadSBarry Smith /* Allocate the column indices */ 884b2f1ab58SBarry Smith retval.row_nnz[i] = nnz; 885b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&retval.icols[i]);CHKERRQ(ierr); 886b2f1ab58SBarry Smith 887c328eeadSBarry Smith /* Administrate the column indices */ 888b2f1ab58SBarry Smith inz = 0; 8896322f4bdSBarry Smith for (j=i; j<maxmrk; j++) { 8906322f4bdSBarry Smith if (iwork[j]) { 891b2f1ab58SBarry Smith retval.icols[i][inz] = j-i; 892b2f1ab58SBarry Smith inz++; 893b2f1ab58SBarry Smith iwork[j]=0; 894b2f1ab58SBarry Smith } 895b2f1ab58SBarry Smith } 896b2f1ab58SBarry Smith retval.nnz += nnz; 897b2f1ab58SBarry Smith }; 898b2f1ab58SBarry Smith ierr = PetscFree(iwork);CHKERRQ(ierr); 899b2f1ab58SBarry Smith *result = retval; 900b2f1ab58SBarry Smith PetscFunctionReturn(0); 901b2f1ab58SBarry Smith } 902b2f1ab58SBarry Smith 903b2f1ab58SBarry Smith 904b2f1ab58SBarry Smith 905b2f1ab58SBarry Smith /* 906b2f1ab58SBarry Smith spbas_keep_upper: 907b2f1ab58SBarry Smith remove the lower part of the matrix: keep the upper part 908b2f1ab58SBarry Smith */ 909b2f1ab58SBarry Smith #undef __FUNCT__ 910b2f1ab58SBarry Smith #define __FUNCT__ "spbas_keep_upper" 911b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper( spbas_matrix * inout_matrix) 912b2f1ab58SBarry Smith { 913b2f1ab58SBarry Smith PetscInt i, j; 914b2f1ab58SBarry Smith PetscInt jstart; 915b2f1ab58SBarry Smith 916b2f1ab58SBarry Smith PetscFunctionBegin; 9176322f4bdSBarry Smith if (inout_matrix->block_data) SETERRQ( PETSC_ERR_SUP_SYS, "Not yet for block data matrices\n"); 9186322f4bdSBarry Smith for (i=0; i<inout_matrix->nrows; i++) { 9196322f4bdSBarry Smith for (jstart=0; (jstart<inout_matrix->row_nnz[i]) && (inout_matrix->icols[i][jstart]<0); jstart++){} 9206322f4bdSBarry Smith if (jstart>0) { 9216322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9226322f4bdSBarry Smith inout_matrix->icols[i][j] = inout_matrix->icols[i][j+jstart]; 923b2f1ab58SBarry Smith } 924b2f1ab58SBarry Smith 9256322f4bdSBarry Smith if (inout_matrix->values) { 9266322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9276322f4bdSBarry Smith inout_matrix->values[i][j] = inout_matrix->values[i][j+jstart]; 928b2f1ab58SBarry Smith } 929b2f1ab58SBarry Smith } 930b2f1ab58SBarry Smith 931b2f1ab58SBarry Smith inout_matrix->row_nnz[i] -= jstart; 932b2f1ab58SBarry Smith 9336322f4bdSBarry Smith inout_matrix->icols[i] = (PetscInt*) realloc((void*) inout_matrix->icols[i], inout_matrix->row_nnz[i]*sizeof(PetscInt)); 934b2f1ab58SBarry Smith 9356322f4bdSBarry Smith if (inout_matrix->values) { 9366322f4bdSBarry Smith inout_matrix->values[i] = (PetscScalar*) realloc((void*) inout_matrix->values[i], inout_matrix->row_nnz[i]*sizeof(PetscScalar)); 937b2f1ab58SBarry Smith } 938b2f1ab58SBarry Smith inout_matrix->nnz -= jstart; 939b2f1ab58SBarry Smith } 940b2f1ab58SBarry Smith } 941b2f1ab58SBarry Smith PetscFunctionReturn(0); 942b2f1ab58SBarry Smith } 943b2f1ab58SBarry Smith 944