1c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas.h> 3b2f1ab58SBarry Smith 4845006b9SBarry Smith /*MC 52692d6eeSBarry Smith MATSOLVERBAS - Provides ICC(k) with drop tolerance 6845006b9SBarry Smith 7845006b9SBarry Smith Works with MATAIJ matrices 8845006b9SBarry Smith 9845006b9SBarry Smith Options Database Keys: 10845006b9SBarry Smith + -pc_factor_levels <l> 11b7c853c4SBarry Smith - -pc_factor_drop_tolerance 12845006b9SBarry Smith 13845006b9SBarry Smith Level: beginner 14845006b9SBarry Smith 15845006b9SBarry Smith Contributed by: Bas van 't Hof 16845006b9SBarry Smith 17b7c853c4SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, PCFactorSetLevels(), PCFactorSetDropTolerance() 18845006b9SBarry Smith 19845006b9SBarry Smith M*/ 209ccc4a9bSBarry Smith 21b2f1ab58SBarry Smith /* 22b2f1ab58SBarry Smith spbas_memory_requirement: 23b2f1ab58SBarry Smith Calculate the number of bytes needed to store tha matrix 24b2f1ab58SBarry Smith */ 25b2f1ab58SBarry Smith #undef __FUNCT__ 26b2f1ab58SBarry Smith #define __FUNCT__ "spbas_memory_requirement" 27b2f1ab58SBarry Smith long int spbas_memory_requirement(spbas_matrix matrix) 28b2f1ab58SBarry Smith { 29c328eeadSBarry Smith long int memreq = 6 * sizeof(PetscInt) + /* nrows, ncols, nnz, n_alloc_icol, 30c328eeadSBarry Smith n_alloc_val, col_idx_type */ 31ace3abfcSBarry Smith sizeof(PetscBool) + /* 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" 61ace3abfcSBarry Smith PetscErrorCode spbas_allocate_pattern(spbas_matrix * result, PetscBool do_values) 62b2f1ab58SBarry Smith { 63b2f1ab58SBarry Smith PetscErrorCode ierr; 64b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 65b2f1ab58SBarry Smith PetscInt col_idx_type = result->col_idx_type; 66b2f1ab58SBarry Smith 67b2f1ab58SBarry Smith PetscFunctionBegin; 68c328eeadSBarry Smith /* Allocate sparseness pattern */ 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; 105ace3abfcSBarry Smith PetscBool do_values = (result->values != PETSC_NULL) ? PETSC_TRUE : PETSC_FALSE; 106ace3abfcSBarry Smith PetscBool 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 ierr = PetscMalloc(nrows*sizeof(PetscInt),&ialloc);CHKERRQ(ierr); 201b2f1ab58SBarry Smith 202b2f1ab58SBarry Smith ihlp1 = ialloc; 203b2f1ab58SBarry Smith ihlp2 = isort; 204b2f1ab58SBarry Smith 205c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 2064e1ff37aSBarry Smith for (istep=1; istep<nrows; istep*=2) { 207c328eeadSBarry Smith /* 208c328eeadSBarry Smith Combine sorted parts 209c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 210c328eeadSBarry Smith of ihlp2 and vhlp2 211c328eeadSBarry Smith 212c328eeadSBarry Smith into one sorted part 213c328eeadSBarry Smith istart:istart+2*istep-1 214c328eeadSBarry Smith of ihlp1 and vhlp1 215c328eeadSBarry Smith */ 2164e1ff37aSBarry Smith for (istart=0; istart<nrows; istart+=2*istep) { 217c328eeadSBarry Smith /* Set counters and bound array part endings */ 218b2f1ab58SBarry Smith i1=istart; i1end = i1+istep; if (i1end>nrows) {i1end=nrows;} 219b2f1ab58SBarry Smith i2=istart+istep; i2end = i2+istep; if (i2end>nrows) {i2end=nrows;} 220b2f1ab58SBarry Smith 221c328eeadSBarry Smith /* Merge the two array parts */ 2224e1ff37aSBarry Smith for (i=istart; i<i2end; i++) { 2234e1ff37aSBarry Smith if (i1<i1end && i2<i2end && spbas_row_order_icol(ihlp2[i1], ihlp2[i2], irow_in, icol_in, col_idx_type) < 0) { 224b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 225b2f1ab58SBarry Smith i1++; 2264e1ff37aSBarry Smith } else if (i2<i2end) { 227b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 228b2f1ab58SBarry Smith i2++; 2294e1ff37aSBarry Smith } else { 230b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 231b2f1ab58SBarry Smith i1++; 232b2f1ab58SBarry Smith } 233b2f1ab58SBarry Smith } 234b2f1ab58SBarry Smith } 235b2f1ab58SBarry Smith 236c328eeadSBarry Smith /* Swap the two array sets */ 237b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 238b2f1ab58SBarry Smith } 239b2f1ab58SBarry Smith 240c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 2414e1ff37aSBarry Smith if (ihlp2 != isort) { 242b2f1ab58SBarry Smith for (i=0; i<nrows; i++) { isort[i] = ihlp2[i]; } 243b2f1ab58SBarry Smith } 244b2f1ab58SBarry Smith ierr = PetscFree(ialloc);CHKERRQ(ierr); 245b2f1ab58SBarry Smith PetscFunctionReturn(0); 246b2f1ab58SBarry Smith } 247b2f1ab58SBarry Smith 248b2f1ab58SBarry Smith 249b2f1ab58SBarry Smith 250b2f1ab58SBarry Smith /* 251b2f1ab58SBarry Smith spbas_compress_pattern: 252b2f1ab58SBarry Smith calculate a compressed sparseness pattern for a sparseness pattern 253b2f1ab58SBarry Smith given in compressed row storage. The compressed sparseness pattern may 254b2f1ab58SBarry Smith require (much) less memory. 255b2f1ab58SBarry Smith */ 256b2f1ab58SBarry Smith #undef __FUNCT__ 257b2f1ab58SBarry Smith #define __FUNCT__ "spbas_compress_pattern" 2584e1ff37aSBarry 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) 259b2f1ab58SBarry Smith { 260b2f1ab58SBarry Smith PetscInt nnz = irow_in[nrows]; 261b2f1ab58SBarry Smith long int mem_orig = (nrows + nnz) * sizeof(PetscInt); 262b2f1ab58SBarry Smith long int mem_compressed; 263b2f1ab58SBarry Smith PetscErrorCode ierr; 264b2f1ab58SBarry Smith PetscInt *isort; 265b2f1ab58SBarry Smith PetscInt *icols; 266b2f1ab58SBarry Smith PetscInt row_nnz; 267b2f1ab58SBarry Smith PetscInt *ipoint; 268ace3abfcSBarry Smith PetscBool *used; 269b2f1ab58SBarry Smith PetscInt ptr; 270b2f1ab58SBarry Smith PetscInt i,j; 271ace3abfcSBarry Smith const PetscBool no_values = PETSC_FALSE; 272b2f1ab58SBarry Smith 273b2f1ab58SBarry Smith PetscFunctionBegin; 274c328eeadSBarry Smith /* Allocate the structure of the new matrix */ 275b2f1ab58SBarry Smith B->nrows = nrows; 276b2f1ab58SBarry Smith B->ncols = ncols; 277b2f1ab58SBarry Smith B->nnz = nnz; 278b2f1ab58SBarry Smith B->col_idx_type= col_idx_type; 279b2f1ab58SBarry Smith B->block_data = PETSC_TRUE; 280b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(B, no_values);CHKERRQ(ierr); 281b2f1ab58SBarry Smith 282c328eeadSBarry Smith /* When using an offset array, set it */ 2834e1ff37aSBarry Smith if (col_idx_type==SPBAS_OFFSET_ARRAY) { 284b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {B->icol0[i] = icol_in[irow_in[i]];} 285b2f1ab58SBarry Smith } 286b2f1ab58SBarry Smith 287c328eeadSBarry Smith /* Allocate the ordering for the rows */ 288b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&isort);CHKERRQ(ierr); 289b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&ipoint);CHKERRQ(ierr); 290ace3abfcSBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscBool),&used);CHKERRQ(ierr); 291b2f1ab58SBarry Smith 292c328eeadSBarry Smith /* Initialize the sorting */ 293ace3abfcSBarry Smith ierr = PetscMemzero((void*) used, nrows*sizeof(PetscBool));CHKERRQ(ierr); 2944e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 295b2f1ab58SBarry Smith B->row_nnz[i] = irow_in[i+1]-irow_in[i]; 296b2f1ab58SBarry Smith isort[i] = i; 297b2f1ab58SBarry Smith ipoint[i]= i; 298b2f1ab58SBarry Smith } 299b2f1ab58SBarry Smith 300c328eeadSBarry Smith /* Sort the rows so that identical columns will be next to each other */ 301b2f1ab58SBarry Smith ierr = spbas_mergesort_icols(nrows, irow_in, icol_in, col_idx_type, isort);CHKERRQ(ierr); 30271d55d18SBarry Smith ierr = PetscInfo(PETSC_NULL,"Rows have been sorted for patterns\n");CHKERRQ(ierr); 303b2f1ab58SBarry Smith 304c328eeadSBarry Smith /* Replace identical rows with the first one in the list */ 3054e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 3064e1ff37aSBarry Smith if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) { 307b2f1ab58SBarry Smith ipoint[isort[i]] = ipoint[isort[i-1]]; 308b2f1ab58SBarry Smith } 309b2f1ab58SBarry Smith } 310b2f1ab58SBarry Smith 311c328eeadSBarry Smith /* Collect the rows which are used*/ 312b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {used[ipoint[i]] = PETSC_TRUE;} 313b2f1ab58SBarry Smith 314c328eeadSBarry Smith /* Calculate needed memory */ 315b2f1ab58SBarry Smith B->n_alloc_icol = 0; 3164e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 317b2f1ab58SBarry Smith if (used[i]) {B->n_alloc_icol += B->row_nnz[i];} 318b2f1ab58SBarry Smith } 31971d55d18SBarry Smith ierr = PetscMalloc(B->n_alloc_icol*sizeof(PetscInt),&B->alloc_icol);CHKERRQ(ierr); 320b2f1ab58SBarry Smith 321c328eeadSBarry Smith /* Fill in the diagonal offsets for the rows which store their own data */ 322b2f1ab58SBarry Smith ptr = 0; 3234e1ff37aSBarry Smith for (i=0; i<B->nrows; i++) { 3244e1ff37aSBarry Smith if (used[i]) { 325b2f1ab58SBarry Smith B->icols[i] = &B->alloc_icol[ptr]; 326b2f1ab58SBarry Smith icols = &icol_in[irow_in[i]]; 327b2f1ab58SBarry Smith row_nnz = B->row_nnz[i]; 3284e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 3294e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 330b2f1ab58SBarry Smith B->icols[i][j] = icols[j]; 331b2f1ab58SBarry Smith } 3324e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 3334e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 334b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-i; 335b2f1ab58SBarry Smith } 3364e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 3374e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 338b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-icols[0]; 339b2f1ab58SBarry Smith } 340b2f1ab58SBarry Smith } 341b2f1ab58SBarry Smith ptr += B->row_nnz[i]; 342b2f1ab58SBarry Smith } 343b2f1ab58SBarry Smith } 344b2f1ab58SBarry Smith 345c328eeadSBarry Smith /* Point to the right places for all data */ 3464e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 347b2f1ab58SBarry Smith B->icols[i] = B->icols[ipoint[i]]; 348b2f1ab58SBarry Smith } 349c328eeadSBarry Smith ierr = PetscInfo(PETSC_NULL,"Row patterns have been compressed\n");CHKERRQ(ierr); 350c328eeadSBarry Smith ierr = PetscInfo1(PETSC_NULL," (%G nonzeros per row)\n", (PetscReal) nnz / (PetscReal) nrows);CHKERRQ(ierr); 351b2f1ab58SBarry Smith 352b2f1ab58SBarry Smith ierr=PetscFree(isort);CHKERRQ(ierr); 353b2f1ab58SBarry Smith ierr=PetscFree(used);CHKERRQ(ierr); 354b2f1ab58SBarry Smith ierr=PetscFree(ipoint);CHKERRQ(ierr); 355b2f1ab58SBarry Smith 356b2f1ab58SBarry Smith mem_compressed = spbas_memory_requirement(*B); 3574e1ff37aSBarry Smith *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig; 358b2f1ab58SBarry Smith PetscFunctionReturn(0); 359b2f1ab58SBarry Smith } 360b2f1ab58SBarry Smith 361b2f1ab58SBarry Smith /* 362b2f1ab58SBarry Smith spbas_incomplete_cholesky 363b2f1ab58SBarry Smith Incomplete Cholesky decomposition 364b2f1ab58SBarry Smith */ 365c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas_cholesky.h> 366b2f1ab58SBarry Smith 367b2f1ab58SBarry Smith /* 368b2f1ab58SBarry Smith spbas_delete : de-allocate the arrays owned by this matrix 369b2f1ab58SBarry Smith */ 370b2f1ab58SBarry Smith #undef __FUNCT__ 371b2f1ab58SBarry Smith #define __FUNCT__ "spbas_delete" 372b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix matrix) 373b2f1ab58SBarry Smith { 374b2f1ab58SBarry Smith PetscInt i; 375b2f1ab58SBarry Smith PetscErrorCode ierr; 3769ccc4a9bSBarry Smith 377b2f1ab58SBarry Smith PetscFunctionBegin; 3789ccc4a9bSBarry Smith if (matrix.block_data) { 379cb9801acSJed Brown ierr=PetscFree(matrix.alloc_icol);CHKERRQ(ierr); 380b2f1ab58SBarry Smith if (matrix.values) {ierr=PetscFree(matrix.alloc_val);CHKERRQ(ierr);} 3819ccc4a9bSBarry Smith } else { 3829ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.icols[i]);CHKERRQ(ierr);} 383b2f1ab58SBarry Smith ierr = PetscFree(matrix.icols);CHKERRQ(ierr); 3849ccc4a9bSBarry Smith if (matrix.values) { 3859ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.values[i]);CHKERRQ(ierr);} 386b2f1ab58SBarry Smith } 387b2f1ab58SBarry Smith } 388b2f1ab58SBarry Smith 389b2f1ab58SBarry Smith ierr=PetscFree(matrix.row_nnz);CHKERRQ(ierr); 390b2f1ab58SBarry Smith ierr=PetscFree(matrix.icols);CHKERRQ(ierr); 3919ccc4a9bSBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) {ierr=PetscFree(matrix.icol0);CHKERRQ(ierr);} 392c31cb41cSBarry Smith ierr=PetscFree(matrix.values);CHKERRQ(ierr); 393b2f1ab58SBarry Smith PetscFunctionReturn(0); 394b2f1ab58SBarry Smith } 395b2f1ab58SBarry Smith 396b2f1ab58SBarry Smith /* 397b2f1ab58SBarry Smith spbas_matrix_to_crs: 398b2f1ab58SBarry Smith Convert an spbas_matrix to compessed row storage 399b2f1ab58SBarry Smith */ 400b2f1ab58SBarry Smith #undef __FUNCT__ 401b2f1ab58SBarry Smith #define __FUNCT__ "spbas_matrix_to_crs" 402b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out) 403b2f1ab58SBarry Smith { 404b2f1ab58SBarry Smith PetscInt nrows = matrix_A.nrows; 405b2f1ab58SBarry Smith PetscInt nnz = matrix_A.nnz; 406b2f1ab58SBarry Smith PetscInt i,j,r_nnz,i0; 407b2f1ab58SBarry Smith PetscInt *irow; 408b2f1ab58SBarry Smith PetscInt *icol; 409b2f1ab58SBarry Smith PetscInt *icol_A; 410b2f1ab58SBarry Smith MatScalar *val; 411b2f1ab58SBarry Smith PetscScalar *val_A; 412b2f1ab58SBarry Smith PetscInt col_idx_type = matrix_A.col_idx_type; 413ace3abfcSBarry Smith PetscBool do_values = matrix_A.values ? PETSC_TRUE : PETSC_FALSE; 414b2f1ab58SBarry Smith PetscErrorCode ierr; 415b2f1ab58SBarry Smith 416b2f1ab58SBarry Smith PetscFunctionBegin; 417b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscInt) * (nrows+1), &irow);CHKERRQ(ierr); 418b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscInt) * nnz, &icol);CHKERRQ(ierr); 4199ccc4a9bSBarry Smith *icol_out = icol; 4209ccc4a9bSBarry Smith *irow_out=irow; 4219ccc4a9bSBarry Smith if (do_values) { 422b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(MatScalar) * nnz, &val);CHKERRQ(ierr); 423b2f1ab58SBarry Smith *val_out = val; *icol_out = icol; *irow_out=irow; 424b2f1ab58SBarry Smith } 425b2f1ab58SBarry Smith 426b2f1ab58SBarry Smith irow[0]=0; 4279ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 428b2f1ab58SBarry Smith r_nnz = matrix_A.row_nnz[i]; 429b2f1ab58SBarry Smith i0 = irow[i]; 430b2f1ab58SBarry Smith irow[i+1] = i0 + r_nnz; 431b2f1ab58SBarry Smith icol_A = matrix_A.icols[i]; 432b2f1ab58SBarry Smith 4339ccc4a9bSBarry Smith if (do_values) { 434b2f1ab58SBarry Smith val_A = matrix_A.values[i]; 4359ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 436b2f1ab58SBarry Smith icol[i0+j] = icol_A[j]; 437b2f1ab58SBarry Smith val[i0+j] = val_A[j]; 438b2f1ab58SBarry Smith } 4399ccc4a9bSBarry Smith } else { 440b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] = icol_A[j]; } 441b2f1ab58SBarry Smith } 442b2f1ab58SBarry Smith 4439ccc4a9bSBarry Smith if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 444b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] += i; } 445b2f1ab58SBarry Smith } 4469ccc4a9bSBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 447b2f1ab58SBarry Smith i0 = matrix_A.icol0[i]; 448b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] += i0; } 449b2f1ab58SBarry Smith } 450b2f1ab58SBarry Smith } 451b2f1ab58SBarry Smith PetscFunctionReturn(0); 452b2f1ab58SBarry Smith } 453b2f1ab58SBarry Smith 454b2f1ab58SBarry Smith 455b2f1ab58SBarry Smith /* 456b2f1ab58SBarry Smith spbas_transpose 457b2f1ab58SBarry Smith return the transpose of a matrix 458b2f1ab58SBarry Smith */ 459b2f1ab58SBarry Smith #undef __FUNCT__ 460b2f1ab58SBarry Smith #define __FUNCT__ "spbas_transpose" 461b2f1ab58SBarry Smith PetscErrorCode spbas_transpose(spbas_matrix in_matrix, spbas_matrix * result) 462b2f1ab58SBarry Smith { 463b2f1ab58SBarry Smith PetscInt col_idx_type = in_matrix.col_idx_type; 464b2f1ab58SBarry Smith PetscInt nnz = in_matrix.nnz; 465b2f1ab58SBarry Smith PetscInt ncols = in_matrix.nrows; 466b2f1ab58SBarry Smith PetscInt nrows = in_matrix.ncols; 467b2f1ab58SBarry Smith PetscInt i,j,k; 468b2f1ab58SBarry Smith PetscInt r_nnz; 469b2f1ab58SBarry Smith PetscInt *irow; 4704efc9174SBarry Smith PetscInt icol0 = 0; 471b2f1ab58SBarry Smith PetscScalar * val; 472b2f1ab58SBarry Smith PetscErrorCode ierr; 4734e1ff37aSBarry Smith 474b2f1ab58SBarry Smith PetscFunctionBegin; 475c328eeadSBarry Smith /* Copy input values */ 476b2f1ab58SBarry Smith result->nrows = nrows; 477b2f1ab58SBarry Smith result->ncols = ncols; 478b2f1ab58SBarry Smith result->nnz = nnz; 479b2f1ab58SBarry Smith result->col_idx_type = SPBAS_COLUMN_NUMBERS; 480b2f1ab58SBarry Smith result->block_data = PETSC_TRUE; 481b2f1ab58SBarry Smith 482c328eeadSBarry Smith /* Allocate sparseness pattern */ 48371d55d18SBarry Smith ierr = spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 484b2f1ab58SBarry Smith 485c328eeadSBarry Smith /* Count the number of nonzeros in each row */ 486b2f1ab58SBarry Smith for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; } 487b2f1ab58SBarry Smith 4889ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 489b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 490b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 4919ccc4a9bSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 492b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[irow[j]]++; } 4939ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 494b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[i+irow[j]]++; } 4959ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 496b2f1ab58SBarry Smith icol0=in_matrix.icol0[i]; 497b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[icol0+irow[j]]++; } 498b2f1ab58SBarry Smith } 499b2f1ab58SBarry Smith } 500b2f1ab58SBarry Smith 501c328eeadSBarry Smith /* Set the pointers to the data */ 502b2f1ab58SBarry Smith ierr = spbas_allocate_data(result);CHKERRQ(ierr); 503b2f1ab58SBarry Smith 504c328eeadSBarry Smith /* Reset the number of nonzeros in each row */ 505b2f1ab58SBarry Smith for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; } 506b2f1ab58SBarry Smith 507c328eeadSBarry Smith /* Fill the data arrays */ 5089ccc4a9bSBarry Smith if (in_matrix.values) { 5099ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 510b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 511b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 512b2f1ab58SBarry Smith val = in_matrix.values[i]; 513b2f1ab58SBarry Smith 514b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) {icol0=0;} 515b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;} 5164efc9174SBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) {icol0=in_matrix.icol0[i];} 5179ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 518b2f1ab58SBarry Smith k = icol0 + irow[j]; 519b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 520b2f1ab58SBarry Smith result->values[k][result->row_nnz[k]] = val[j]; 521b2f1ab58SBarry Smith result->row_nnz[k]++; 522b2f1ab58SBarry Smith } 523b2f1ab58SBarry Smith } 5249ccc4a9bSBarry Smith } else { 5259ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 526b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 527b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 528b2f1ab58SBarry Smith 529b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) {icol0=0;} 530b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;} 5319ccc4a9bSBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) {icol0=in_matrix.icol0[i];} 532b2f1ab58SBarry Smith 5339ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 534b2f1ab58SBarry Smith k = icol0 + irow[j]; 535b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 536b2f1ab58SBarry Smith result->row_nnz[k]++; 537b2f1ab58SBarry Smith } 538b2f1ab58SBarry Smith } 539b2f1ab58SBarry Smith } 540b2f1ab58SBarry Smith PetscFunctionReturn(0); 541b2f1ab58SBarry Smith } 542b2f1ab58SBarry Smith 543b2f1ab58SBarry Smith /* 544b2f1ab58SBarry Smith spbas_mergesort 545b2f1ab58SBarry Smith 546b2f1ab58SBarry Smith mergesort for an array of intergers and an array of associated 547b2f1ab58SBarry Smith reals 548b2f1ab58SBarry Smith 549b2f1ab58SBarry Smith on output, icol[0..nnz-1] is increasing; 550b2f1ab58SBarry Smith val[0..nnz-1] has undergone the same permutation as icol 551b2f1ab58SBarry Smith 552c328eeadSBarry Smith NB: val may be PETSC_NULL: in that case, only the integers are sorted 553b2f1ab58SBarry Smith 554b2f1ab58SBarry Smith */ 555b2f1ab58SBarry Smith #undef __FUNCT__ 556b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort" 557b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val) 558b2f1ab58SBarry Smith { 559c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 560c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 561c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 562c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 563c328eeadSBarry Smith PetscScalar *valloc=PETSC_NULL; 564c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 565b2f1ab58SBarry Smith PetscScalar *vswap; 566c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 567c328eeadSBarry Smith PetscScalar *vhlp1=PETSC_NULL; /* (arrays under construction) */ 568c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 569c328eeadSBarry Smith PetscScalar *vhlp2=PETSC_NULL; 570b2f1ab58SBarry Smith PetscErrorCode ierr; 571b2f1ab58SBarry Smith 572b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&ialloc);CHKERRQ(ierr); 573b2f1ab58SBarry Smith ihlp1 = ialloc; 574b2f1ab58SBarry Smith ihlp2 = icol; 575b2f1ab58SBarry Smith 5769ccc4a9bSBarry Smith if (val) { 577b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscScalar),&valloc);CHKERRQ(ierr); 578b2f1ab58SBarry Smith vhlp1 = valloc; 579b2f1ab58SBarry Smith vhlp2 = val; 580b2f1ab58SBarry Smith } 581b2f1ab58SBarry Smith 582b2f1ab58SBarry Smith 583c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 5849ccc4a9bSBarry Smith for (istep=1; istep<nnz; istep*=2) { 585c328eeadSBarry Smith /* 586c328eeadSBarry Smith Combine sorted parts 587c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 588c328eeadSBarry Smith of ihlp2 and vhlp2 589c328eeadSBarry Smith 590c328eeadSBarry Smith into one sorted part 591c328eeadSBarry Smith istart:istart+2*istep-1 592c328eeadSBarry Smith of ihlp1 and vhlp1 593c328eeadSBarry Smith */ 5949ccc4a9bSBarry Smith for (istart=0; istart<nnz; istart+=2*istep) { 595c328eeadSBarry Smith /* Set counters and bound array part endings */ 596b2f1ab58SBarry Smith i1=istart; i1end = i1+istep; if (i1end>nnz) {i1end=nnz;} 597b2f1ab58SBarry Smith i2=istart+istep; i2end = i2+istep; if (i2end>nnz) {i2end=nnz;} 598b2f1ab58SBarry Smith 599c328eeadSBarry Smith /* Merge the two array parts */ 6009ccc4a9bSBarry Smith if (val) { 6019ccc4a9bSBarry Smith for (i=istart; i<i2end; i++) { 6029ccc4a9bSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 603b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 604b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 605b2f1ab58SBarry Smith i1++; 6069ccc4a9bSBarry Smith } else if (i2<i2end) { 607b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 608b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i2]; 609b2f1ab58SBarry Smith i2++; 6109ccc4a9bSBarry Smith } else { 611b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 612b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 613b2f1ab58SBarry Smith i1++; 614b2f1ab58SBarry Smith } 615b2f1ab58SBarry Smith } 6169ccc4a9bSBarry Smith } else { 6176322f4bdSBarry Smith for (i=istart; i<i2end; i++) { 6186322f4bdSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 619b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 620b2f1ab58SBarry Smith i1++; 6216322f4bdSBarry Smith } else if (i2<i2end) { 622b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 623b2f1ab58SBarry Smith i2++; 6246322f4bdSBarry Smith } else { 625b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 626b2f1ab58SBarry Smith i1++; 627b2f1ab58SBarry Smith } 628b2f1ab58SBarry Smith } 629b2f1ab58SBarry Smith } 630b2f1ab58SBarry Smith } 631b2f1ab58SBarry Smith 632c328eeadSBarry Smith /* Swap the two array sets */ 633b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 634b2f1ab58SBarry Smith vswap = vhlp2; vhlp2 = vhlp1; vhlp1 = vswap; 635b2f1ab58SBarry Smith } 636b2f1ab58SBarry Smith 637c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 6386322f4bdSBarry Smith if (ihlp2 != icol) { 639b2f1ab58SBarry Smith for (i=0; i<nnz; i++) { icol[i] = ihlp2[i]; } 6406322f4bdSBarry Smith if (val) { 641b2f1ab58SBarry Smith for (i=0; i<nnz; i++) { val[i] = vhlp2[i]; } 642b2f1ab58SBarry Smith } 643b2f1ab58SBarry Smith } 644b2f1ab58SBarry Smith 645b2f1ab58SBarry Smith ierr = PetscFree(ialloc);CHKERRQ(ierr); 646b2f1ab58SBarry Smith if (val) {ierr = PetscFree(valloc);CHKERRQ(ierr);} 647b2f1ab58SBarry Smith PetscFunctionReturn(0); 648b2f1ab58SBarry Smith } 649b2f1ab58SBarry Smith 650b2f1ab58SBarry Smith /* 651b2f1ab58SBarry Smith spbas_apply_reordering_rows: 652b2f1ab58SBarry Smith apply the given reordering to the rows: matrix_A = matrix_A(perm,:); 653b2f1ab58SBarry Smith */ 654b2f1ab58SBarry Smith #undef __FUNCT__ 655b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_rows" 656b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_rows(spbas_matrix *matrix_A, const PetscInt *permutation) 657b2f1ab58SBarry Smith { 658b2f1ab58SBarry Smith PetscInt i,j,ip; 659b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 660b2f1ab58SBarry Smith PetscInt * row_nnz; 661b2f1ab58SBarry Smith PetscInt **icols; 662ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 663c328eeadSBarry Smith PetscScalar **vals=PETSC_NULL; 664b2f1ab58SBarry Smith PetscErrorCode ierr; 665b2f1ab58SBarry Smith 666b2f1ab58SBarry Smith PetscFunctionBegin; 667e32f2f54SBarry Smith if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern\n"); 668b2f1ab58SBarry Smith 6696322f4bdSBarry Smith if (do_values) { 670b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscScalar*)*nrows, &vals);CHKERRQ(ierr); 671b2f1ab58SBarry Smith } 672b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*nrows, &row_nnz);CHKERRQ(ierr); 673b2f1ab58SBarry Smith ierr = PetscMalloc(sizeof(PetscInt*)*nrows, &icols);CHKERRQ(ierr); 674b2f1ab58SBarry Smith 6756322f4bdSBarry Smith for (i=0; i<nrows;i++) { 676b2f1ab58SBarry Smith ip = permutation[i]; 677b2f1ab58SBarry Smith if (do_values) {vals[i] = matrix_A->values[ip];} 678b2f1ab58SBarry Smith icols[i] = matrix_A->icols[ip]; 679b2f1ab58SBarry Smith row_nnz[i] = matrix_A->row_nnz[ip]; 680b2f1ab58SBarry Smith for (j=0; j<row_nnz[i]; j++) { icols[i][j] += ip-i; } 681b2f1ab58SBarry Smith } 682b2f1ab58SBarry Smith 683b2f1ab58SBarry Smith if (do_values) { ierr = PetscFree(matrix_A->values);CHKERRQ(ierr);} 684b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->icols);CHKERRQ(ierr); 685b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->row_nnz);CHKERRQ(ierr); 686b2f1ab58SBarry Smith 687b2f1ab58SBarry Smith if (do_values) { matrix_A->values = vals; } 688b2f1ab58SBarry Smith matrix_A->icols = icols; 689b2f1ab58SBarry Smith matrix_A->row_nnz = row_nnz; 690b2f1ab58SBarry Smith 691b2f1ab58SBarry Smith PetscFunctionReturn(0); 692b2f1ab58SBarry Smith } 693b2f1ab58SBarry Smith 694b2f1ab58SBarry Smith 695b2f1ab58SBarry Smith /* 696b2f1ab58SBarry Smith spbas_apply_reordering_cols: 697b2f1ab58SBarry Smith apply the given reordering to the columns: matrix_A(:,perm) = matrix_A; 698b2f1ab58SBarry Smith */ 699b2f1ab58SBarry Smith #undef __FUNCT__ 700b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_cols" 701b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_cols(spbas_matrix *matrix_A,const PetscInt *permutation) 702b2f1ab58SBarry Smith { 703b2f1ab58SBarry Smith PetscInt i,j; 704b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 705b2f1ab58SBarry Smith PetscInt row_nnz; 706b2f1ab58SBarry Smith PetscInt *icols; 707ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 708c328eeadSBarry Smith PetscScalar *vals=PETSC_NULL; 709b2f1ab58SBarry Smith PetscErrorCode ierr; 710b2f1ab58SBarry Smith 711b2f1ab58SBarry Smith PetscFunctionBegin; 712e32f2f54SBarry Smith if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "must have diagonal offsets in pattern\n"); 713b2f1ab58SBarry Smith 7146322f4bdSBarry Smith for (i=0; i<nrows;i++) { 715b2f1ab58SBarry Smith icols = matrix_A->icols[i]; 716b2f1ab58SBarry Smith row_nnz = matrix_A->row_nnz[i]; 717b2f1ab58SBarry Smith if (do_values) { vals = matrix_A->values[i]; } 718b2f1ab58SBarry Smith 7196322f4bdSBarry Smith for (j=0; j<row_nnz; j++) { 720b2f1ab58SBarry Smith icols[j] = permutation[i+icols[j]]-i; 721b2f1ab58SBarry Smith } 722b2f1ab58SBarry Smith ierr = spbas_mergesort(row_nnz, icols, vals);CHKERRQ(ierr); 723b2f1ab58SBarry Smith } 724b2f1ab58SBarry Smith PetscFunctionReturn(0); 725b2f1ab58SBarry Smith } 726b2f1ab58SBarry Smith 727b2f1ab58SBarry Smith /* 728b2f1ab58SBarry Smith spbas_apply_reordering: 729b2f1ab58SBarry Smith apply the given reordering: matrix_A(perm,perm) = matrix_A; 730b2f1ab58SBarry Smith */ 731b2f1ab58SBarry Smith #undef __FUNCT__ 732b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering" 733b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering(spbas_matrix *matrix_A, const PetscInt *permutation, const PetscInt * inv_perm) 734b2f1ab58SBarry Smith { 735b2f1ab58SBarry Smith PetscErrorCode ierr; 736*5fd66863SKarl Rupp 737b2f1ab58SBarry Smith PetscFunctionBegin; 738b2f1ab58SBarry Smith ierr = spbas_apply_reordering_rows(matrix_A, inv_perm);CHKERRQ(ierr); 739b2f1ab58SBarry Smith ierr = spbas_apply_reordering_cols(matrix_A, permutation);CHKERRQ(ierr); 740b2f1ab58SBarry Smith PetscFunctionReturn(0); 741b2f1ab58SBarry Smith } 742b2f1ab58SBarry Smith 743b2f1ab58SBarry Smith #undef __FUNCT__ 744b2f1ab58SBarry Smith #define __FUNCT__ "spbas_pattern_only" 745b2f1ab58SBarry Smith PetscErrorCode spbas_pattern_only(PetscInt nrows, PetscInt ncols, PetscInt *ai, PetscInt *aj, spbas_matrix * result) 746b2f1ab58SBarry Smith { 747b2f1ab58SBarry Smith spbas_matrix retval; 748b2f1ab58SBarry Smith PetscInt i, j, i0, r_nnz; 749b2f1ab58SBarry Smith PetscErrorCode ierr; 750b2f1ab58SBarry Smith 751b2f1ab58SBarry Smith PetscFunctionBegin; 752c328eeadSBarry Smith /* Copy input values */ 753b2f1ab58SBarry Smith retval.nrows = nrows; 754b2f1ab58SBarry Smith retval.ncols = ncols; 755b2f1ab58SBarry Smith retval.nnz = ai[nrows]; 756b2f1ab58SBarry Smith 757b2f1ab58SBarry Smith retval.block_data = PETSC_TRUE; 758b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 759b2f1ab58SBarry Smith 760c328eeadSBarry Smith /* Allocate output matrix */ 761b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(&retval, PETSC_FALSE);CHKERRQ(ierr); 762b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {retval.row_nnz[i] = ai[i+1]-ai[i];} 763b2f1ab58SBarry Smith ierr = spbas_allocate_data(&retval);CHKERRQ(ierr); 764c328eeadSBarry Smith /* Copy the structure */ 7656322f4bdSBarry Smith for (i = 0; i<retval.nrows; i++) { 766b2f1ab58SBarry Smith i0 = ai[i]; 767b2f1ab58SBarry Smith r_nnz = ai[i+1]-i0; 768b2f1ab58SBarry Smith 7696322f4bdSBarry Smith for (j=0; j<r_nnz; j++) { 770b2f1ab58SBarry Smith retval.icols[i][j] = aj[i0+j]-i; 771b2f1ab58SBarry Smith } 772b2f1ab58SBarry Smith } 773b2f1ab58SBarry Smith *result = retval; 774b2f1ab58SBarry Smith PetscFunctionReturn(0); 775b2f1ab58SBarry Smith } 776b2f1ab58SBarry Smith 777b2f1ab58SBarry Smith 778b2f1ab58SBarry Smith /* 779b2f1ab58SBarry Smith spbas_mark_row_power: 780b2f1ab58SBarry Smith Mark the columns in row 'row' which are nonzero in 781b2f1ab58SBarry Smith matrix^2log(marker). 782b2f1ab58SBarry Smith */ 783b2f1ab58SBarry Smith #undef __FUNCT__ 784b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mark_row_power" 785be332245SKarl Rupp PetscErrorCode spbas_mark_row_power(PetscInt *iwork, /* marker-vector */ 786c328eeadSBarry Smith PetscInt row, /* row for which the columns are marked */ 787c328eeadSBarry Smith spbas_matrix * in_matrix, /* matrix for which the power is being calculated */ 788c328eeadSBarry Smith PetscInt marker, /* marker-value: 2^power */ 789c328eeadSBarry Smith PetscInt minmrk, /* lower bound for marked points */ 790c328eeadSBarry Smith PetscInt maxmrk) /* upper bound for marked points */ 791b2f1ab58SBarry Smith { 792b2f1ab58SBarry Smith PetscErrorCode ierr; 793b2f1ab58SBarry Smith PetscInt i,j, nnz; 794b2f1ab58SBarry Smith 795b2f1ab58SBarry Smith PetscFunctionBegin; 796b2f1ab58SBarry Smith nnz = in_matrix->row_nnz[row]; 797b2f1ab58SBarry Smith 798c328eeadSBarry Smith /* For higher powers, call this function recursively */ 7996322f4bdSBarry Smith if (marker>1) { 8006322f4bdSBarry Smith for (i=0; i<nnz;i++) { 801b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 8026322f4bdSBarry Smith if (minmrk<=j && j<maxmrk && iwork[j] < marker) { 80371d55d18SBarry Smith ierr = spbas_mark_row_power(iwork, row + in_matrix->icols[row][i],in_matrix, marker/2,minmrk,maxmrk);CHKERRQ(ierr); 804b2f1ab58SBarry Smith iwork[j] |= marker; 805b2f1ab58SBarry Smith } 806b2f1ab58SBarry Smith } 8076322f4bdSBarry Smith } else { 808c328eeadSBarry Smith /* Mark the columns reached */ 8096322f4bdSBarry Smith for (i=0; i<nnz;i++) { 810b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 811b2f1ab58SBarry Smith if (minmrk<=j && j<maxmrk) { iwork[j] |= 1; } 812b2f1ab58SBarry Smith } 813b2f1ab58SBarry Smith } 814b2f1ab58SBarry Smith PetscFunctionReturn(0); 815b2f1ab58SBarry Smith } 816b2f1ab58SBarry Smith 817b2f1ab58SBarry Smith 818b2f1ab58SBarry Smith /* 819b2f1ab58SBarry Smith spbas_power 820b2f1ab58SBarry Smith Calculate sparseness patterns for incomplete Cholesky decompositions 821b2f1ab58SBarry Smith of a given order: (almost) all nonzeros of the matrix^(order+1) which 822b2f1ab58SBarry Smith are inside the band width are found and stored in the output sparseness 823b2f1ab58SBarry Smith pattern. 824b2f1ab58SBarry Smith */ 825b2f1ab58SBarry Smith #undef __FUNCT__ 826b2f1ab58SBarry Smith #define __FUNCT__ "spbas_power" 827b2f1ab58SBarry Smith PetscErrorCode spbas_power (spbas_matrix in_matrix,PetscInt power, spbas_matrix * result) 828b2f1ab58SBarry Smith { 829b2f1ab58SBarry Smith spbas_matrix retval; 830b2f1ab58SBarry Smith PetscInt nrows = in_matrix.nrows; 831b2f1ab58SBarry Smith PetscInt ncols = in_matrix.ncols; 832b2f1ab58SBarry Smith PetscInt i, j, kend; 833b2f1ab58SBarry Smith PetscInt nnz, inz; 834b2f1ab58SBarry Smith PetscInt *iwork; 835b2f1ab58SBarry Smith PetscInt marker; 836b2f1ab58SBarry Smith PetscInt maxmrk=0; 837b2f1ab58SBarry Smith PetscErrorCode ierr; 838b2f1ab58SBarry Smith 839b2f1ab58SBarry Smith PetscFunctionBegin; 840e32f2f54SBarry Smith if (in_matrix.col_idx_type != SPBAS_DIAGONAL_OFFSETS) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS,"must have diagonal offsets in pattern\n"); 841e32f2f54SBarry Smith if (ncols != nrows) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Dimension error\n"); 842e32f2f54SBarry Smith if (in_matrix.values) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input array must be sparseness pattern (no values)"); 843e32f2f54SBarry Smith if (power<=0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Power must be 1 or up"); 844b2f1ab58SBarry Smith 845c328eeadSBarry Smith /* Copy input values*/ 846b2f1ab58SBarry Smith retval.nrows = ncols; 847b2f1ab58SBarry Smith retval.ncols = nrows; 848b2f1ab58SBarry Smith retval.nnz = 0; 849b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 850b2f1ab58SBarry Smith retval.block_data = PETSC_FALSE; 851b2f1ab58SBarry Smith 852c328eeadSBarry Smith /* Allocate sparseness pattern */ 85371d55d18SBarry Smith ierr = spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 854b2f1ab58SBarry Smith 855c328eeadSBarry Smith /* Allocate marker array */ 856b2f1ab58SBarry Smith ierr = PetscMalloc(nrows * sizeof(PetscInt), &iwork);CHKERRQ(ierr); 857b2f1ab58SBarry Smith 858c328eeadSBarry Smith /* Erase the pattern for this row */ 8594e1ff37aSBarry Smith ierr = PetscMemzero((void *) iwork, retval.nrows*sizeof(PetscInt));CHKERRQ(ierr); 860b2f1ab58SBarry Smith 861c328eeadSBarry Smith /* Calculate marker values */ 862b2f1ab58SBarry Smith marker = 1; for (i=1; i<power; i++) {marker*=2;} 863b2f1ab58SBarry Smith 8646322f4bdSBarry Smith for (i=0; i<nrows; i++) { 865c328eeadSBarry Smith /* Calculate the pattern for each row */ 866b2f1ab58SBarry Smith 867b2f1ab58SBarry Smith nnz = in_matrix.row_nnz[i]; 868b2f1ab58SBarry Smith kend = i+in_matrix.icols[i][nnz-1]; 869b2f1ab58SBarry Smith if (maxmrk<=kend) {maxmrk=kend+1;} 87071d55d18SBarry Smith ierr = spbas_mark_row_power(iwork, i, &in_matrix, marker, i, maxmrk);CHKERRQ(ierr); 871b2f1ab58SBarry Smith 872c328eeadSBarry Smith /* Count the columns*/ 873b2f1ab58SBarry Smith nnz = 0; 874b2f1ab58SBarry Smith for (j=i; j<maxmrk; j++) {nnz+= (iwork[j]!=0);} 875b2f1ab58SBarry Smith 876c328eeadSBarry Smith /* Allocate the column indices */ 877b2f1ab58SBarry Smith retval.row_nnz[i] = nnz; 878b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&retval.icols[i]);CHKERRQ(ierr); 879b2f1ab58SBarry Smith 880c328eeadSBarry Smith /* Administrate the column indices */ 881b2f1ab58SBarry Smith inz = 0; 8826322f4bdSBarry Smith for (j=i; j<maxmrk; j++) { 8836322f4bdSBarry Smith if (iwork[j]) { 884b2f1ab58SBarry Smith retval.icols[i][inz] = j-i; 885b2f1ab58SBarry Smith inz++; 886b2f1ab58SBarry Smith iwork[j]=0; 887b2f1ab58SBarry Smith } 888b2f1ab58SBarry Smith } 889b2f1ab58SBarry Smith retval.nnz += nnz; 890b2f1ab58SBarry Smith }; 891b2f1ab58SBarry Smith ierr = PetscFree(iwork);CHKERRQ(ierr); 892b2f1ab58SBarry Smith *result = retval; 893b2f1ab58SBarry Smith PetscFunctionReturn(0); 894b2f1ab58SBarry Smith } 895b2f1ab58SBarry Smith 896b2f1ab58SBarry Smith 897b2f1ab58SBarry Smith 898b2f1ab58SBarry Smith /* 899b2f1ab58SBarry Smith spbas_keep_upper: 900b2f1ab58SBarry Smith remove the lower part of the matrix: keep the upper part 901b2f1ab58SBarry Smith */ 902b2f1ab58SBarry Smith #undef __FUNCT__ 903b2f1ab58SBarry Smith #define __FUNCT__ "spbas_keep_upper" 904b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper(spbas_matrix * inout_matrix) 905b2f1ab58SBarry Smith { 906b2f1ab58SBarry Smith PetscInt i, j; 907b2f1ab58SBarry Smith PetscInt jstart; 908b2f1ab58SBarry Smith 909b2f1ab58SBarry Smith PetscFunctionBegin; 910e32f2f54SBarry Smith if (inout_matrix->block_data) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Not yet for block data matrices\n"); 9116322f4bdSBarry Smith for (i=0; i<inout_matrix->nrows; i++) { 9126322f4bdSBarry Smith for (jstart=0; (jstart<inout_matrix->row_nnz[i]) && (inout_matrix->icols[i][jstart]<0); jstart++) {} 9136322f4bdSBarry Smith if (jstart>0) { 9146322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9156322f4bdSBarry Smith inout_matrix->icols[i][j] = inout_matrix->icols[i][j+jstart]; 916b2f1ab58SBarry Smith } 917b2f1ab58SBarry Smith 9186322f4bdSBarry Smith if (inout_matrix->values) { 9196322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9206322f4bdSBarry Smith inout_matrix->values[i][j] = inout_matrix->values[i][j+jstart]; 921b2f1ab58SBarry Smith } 922b2f1ab58SBarry Smith } 923b2f1ab58SBarry Smith 924b2f1ab58SBarry Smith inout_matrix->row_nnz[i] -= jstart; 925b2f1ab58SBarry Smith 9266322f4bdSBarry Smith inout_matrix->icols[i] = (PetscInt*) realloc((void*) inout_matrix->icols[i], inout_matrix->row_nnz[i]*sizeof(PetscInt)); 927b2f1ab58SBarry Smith 9286322f4bdSBarry Smith if (inout_matrix->values) { 9296322f4bdSBarry Smith inout_matrix->values[i] = (PetscScalar*) realloc((void*) inout_matrix->values[i], inout_matrix->row_nnz[i]*sizeof(PetscScalar)); 930b2f1ab58SBarry Smith } 931b2f1ab58SBarry Smith inout_matrix->nnz -= jstart; 932b2f1ab58SBarry Smith } 933b2f1ab58SBarry Smith } 934b2f1ab58SBarry Smith PetscFunctionReturn(0); 935b2f1ab58SBarry Smith } 936b2f1ab58SBarry Smith 937