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: 10*0053dfc8SBarry Smith + -pc_factor_levels <l> - number of levels of fill 11*0053dfc8SBarry Smith - -pc_factor_drop_tolerance - is not currently hooked up to do anything 12845006b9SBarry Smith 13*0053dfc8SBarry Smith Level: intermediate 14845006b9SBarry Smith 15845006b9SBarry Smith Contributed by: Bas van 't Hof 16845006b9SBarry Smith 17*0053dfc8SBarry Smith Notes: Since this currently hooked up to use drop tolerance it should produce the same factors and hence convergence as the PETSc ICC, for higher 18*0053dfc8SBarry Smith levels of fill it does not. This needs to be investigated. Unless you are interested in drop tolerance ICC and willing to work through the code 19*0053dfc8SBarry Smith we recommend not using this functionality. 20*0053dfc8SBarry Smith 21b7c853c4SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, PCFactorSetLevels(), PCFactorSetDropTolerance() 22845006b9SBarry Smith 23845006b9SBarry Smith M*/ 249ccc4a9bSBarry Smith 25b2f1ab58SBarry Smith /* 26b2f1ab58SBarry Smith spbas_memory_requirement: 27b2f1ab58SBarry Smith Calculate the number of bytes needed to store tha matrix 28b2f1ab58SBarry Smith */ 29b2f1ab58SBarry Smith #undef __FUNCT__ 30b2f1ab58SBarry Smith #define __FUNCT__ "spbas_memory_requirement" 31b2f1ab58SBarry Smith long int spbas_memory_requirement(spbas_matrix matrix) 32b2f1ab58SBarry Smith { 332205254eSKarl Rupp long int memreq = 6 * sizeof(PetscInt) + /* nrows, ncols, nnz, n_alloc_icol, n_alloc_val, col_idx_type */ 34ace3abfcSBarry Smith sizeof(PetscBool) + /* block_data */ 35c328eeadSBarry Smith sizeof(PetscScalar**) + /* values */ 36c328eeadSBarry Smith sizeof(PetscScalar*) + /* alloc_val */ 37c328eeadSBarry Smith 2 * sizeof(int**) + /* icols, icols0 */ 38c328eeadSBarry Smith 2 * sizeof(PetscInt*) + /* row_nnz, alloc_icol */ 39c328eeadSBarry Smith matrix.nrows * sizeof(PetscInt) + /* row_nnz[*] */ 40c328eeadSBarry Smith matrix.nrows * sizeof(PetscInt*); /* icols[*] */ 41b2f1ab58SBarry Smith 42c328eeadSBarry Smith /* icol0[*] */ 432205254eSKarl Rupp if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) memreq += matrix.nrows * sizeof(PetscInt); 44b2f1ab58SBarry Smith 45c328eeadSBarry Smith /* icols[*][*] */ 462205254eSKarl Rupp if (matrix.block_data) memreq += matrix.n_alloc_icol * sizeof(PetscInt); 472205254eSKarl Rupp else memreq += matrix.nnz * sizeof(PetscInt); 48b2f1ab58SBarry Smith 494e1ff37aSBarry Smith if (matrix.values) { 50c328eeadSBarry Smith memreq += matrix.nrows * sizeof(PetscScalar*); /* values[*] */ 51c328eeadSBarry Smith /* values[*][*] */ 522205254eSKarl Rupp if (matrix.block_data) memreq += matrix.n_alloc_val * sizeof(PetscScalar); 532205254eSKarl Rupp else memreq += matrix.nnz * sizeof(PetscScalar); 54b2f1ab58SBarry Smith } 55b2f1ab58SBarry Smith return memreq; 56b2f1ab58SBarry Smith } 57b2f1ab58SBarry Smith 58b2f1ab58SBarry Smith /* 59b2f1ab58SBarry Smith spbas_allocate_pattern: 60b2f1ab58SBarry Smith allocate the pattern arrays row_nnz, icols and optionally values 61b2f1ab58SBarry Smith */ 62b2f1ab58SBarry Smith #undef __FUNCT__ 63b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_pattern" 64ace3abfcSBarry Smith PetscErrorCode spbas_allocate_pattern(spbas_matrix * result, PetscBool do_values) 65b2f1ab58SBarry Smith { 66b2f1ab58SBarry Smith PetscErrorCode ierr; 67b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 68b2f1ab58SBarry Smith PetscInt col_idx_type = result->col_idx_type; 69b2f1ab58SBarry Smith 70b2f1ab58SBarry Smith PetscFunctionBegin; 71c328eeadSBarry Smith /* Allocate sparseness pattern */ 72785e854fSJed Brown ierr = PetscMalloc1(nrows,&result->row_nnz);CHKERRQ(ierr); 73785e854fSJed Brown ierr = PetscMalloc1(nrows,&result->icols);CHKERRQ(ierr); 74b2f1ab58SBarry Smith 75c328eeadSBarry Smith /* If offsets are given wrt an array, create array */ 764e1ff37aSBarry Smith if (col_idx_type == SPBAS_OFFSET_ARRAY) { 77785e854fSJed Brown ierr = PetscMalloc1(nrows,&result->icol0);CHKERRQ(ierr); 784e1ff37aSBarry Smith } else { 790298fd71SBarry Smith result->icol0 = NULL; 80b2f1ab58SBarry Smith } 81b2f1ab58SBarry Smith 82c328eeadSBarry Smith /* If values are given, allocate values array */ 834e1ff37aSBarry Smith if (do_values) { 84785e854fSJed Brown ierr = PetscMalloc1(nrows,&result->values);CHKERRQ(ierr); 854e1ff37aSBarry Smith } else { 860298fd71SBarry Smith result->values = NULL; 87b2f1ab58SBarry Smith } 88b2f1ab58SBarry Smith PetscFunctionReturn(0); 89b2f1ab58SBarry Smith } 90b2f1ab58SBarry Smith 91b2f1ab58SBarry Smith /* 92b2f1ab58SBarry Smith spbas_allocate_data: 93b2f1ab58SBarry Smith in case of block_data: 94b2f1ab58SBarry Smith Allocate the data arrays alloc_icol and optionally alloc_val, 95b2f1ab58SBarry Smith set appropriate pointers from icols and values; 96b2f1ab58SBarry Smith in case of !block_data: 97b2f1ab58SBarry Smith Allocate the arrays icols[i] and optionally values[i] 98b2f1ab58SBarry Smith */ 99b2f1ab58SBarry Smith #undef __FUNCT__ 100b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_data" 101b2f1ab58SBarry Smith PetscErrorCode spbas_allocate_data(spbas_matrix * result) 102b2f1ab58SBarry Smith { 103b2f1ab58SBarry Smith PetscInt i; 104b2f1ab58SBarry Smith PetscInt nnz = result->nnz; 105b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 106b2f1ab58SBarry Smith PetscInt r_nnz; 107b2f1ab58SBarry Smith PetscErrorCode ierr; 1086c4ed002SBarry Smith PetscBool do_values = (result->values) ? PETSC_TRUE : PETSC_FALSE; 109ace3abfcSBarry Smith PetscBool block_data = result->block_data; 110b2f1ab58SBarry Smith 111b2f1ab58SBarry Smith PetscFunctionBegin; 1124e1ff37aSBarry Smith if (block_data) { 113c328eeadSBarry Smith /* Allocate the column number array and point to it */ 114b2f1ab58SBarry Smith result->n_alloc_icol = nnz; 1152205254eSKarl Rupp 116785e854fSJed Brown ierr = PetscMalloc1(nnz, &result->alloc_icol);CHKERRQ(ierr); 1172205254eSKarl Rupp 118b2f1ab58SBarry Smith result->icols[0] = result->alloc_icol; 1194e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 120b2f1ab58SBarry Smith result->icols[i] = result->icols[i-1] + result->row_nnz[i-1]; 121b2f1ab58SBarry Smith } 122b2f1ab58SBarry Smith 123c328eeadSBarry Smith /* Allocate the value array and point to it */ 1244e1ff37aSBarry Smith if (do_values) { 125b2f1ab58SBarry Smith result->n_alloc_val = nnz; 1262205254eSKarl Rupp 127785e854fSJed Brown ierr = PetscMalloc1(nnz, &result->alloc_val);CHKERRQ(ierr); 1282205254eSKarl Rupp 129b2f1ab58SBarry Smith result->values[0] = result->alloc_val; 1304e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 131b2f1ab58SBarry Smith result->values[i] = result->values[i-1] + result->row_nnz[i-1]; 132b2f1ab58SBarry Smith } 133b2f1ab58SBarry Smith } 1344e1ff37aSBarry Smith } else { 1354e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 136b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 137785e854fSJed Brown ierr = PetscMalloc1(r_nnz, &result->icols[i]);CHKERRQ(ierr); 138b2f1ab58SBarry Smith } 1394e1ff37aSBarry Smith if (do_values) { 1404e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 141b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 142785e854fSJed Brown ierr = PetscMalloc1(r_nnz, &result->values[i]);CHKERRQ(ierr); 143b2f1ab58SBarry Smith } 144b2f1ab58SBarry Smith } 145b2f1ab58SBarry Smith } 146b2f1ab58SBarry Smith PetscFunctionReturn(0); 147b2f1ab58SBarry Smith } 148b2f1ab58SBarry Smith 149b2f1ab58SBarry Smith /* 150b2f1ab58SBarry Smith spbas_row_order_icol 151b2f1ab58SBarry Smith determine if row i1 should come 152b2f1ab58SBarry Smith + before row i2 in the sorted rows (return -1), 153b2f1ab58SBarry Smith + after (return 1) 154b2f1ab58SBarry Smith + is identical (return 0). 155b2f1ab58SBarry Smith */ 156b2f1ab58SBarry Smith #undef __FUNCT__ 157b2f1ab58SBarry Smith #define __FUNCT__ "spbas_row_order_icol" 158b2f1ab58SBarry Smith int spbas_row_order_icol(PetscInt i1, PetscInt i2, PetscInt *irow_in, PetscInt *icol_in,PetscInt col_idx_type) 159b2f1ab58SBarry Smith { 160b2f1ab58SBarry Smith PetscInt j; 161b2f1ab58SBarry Smith PetscInt nnz1 = irow_in[i1+1] - irow_in[i1]; 162b2f1ab58SBarry Smith PetscInt nnz2 = irow_in[i2+1] - irow_in[i2]; 163b2f1ab58SBarry Smith PetscInt * icol1 = &icol_in[irow_in[i1]]; 164b2f1ab58SBarry Smith PetscInt * icol2 = &icol_in[irow_in[i2]]; 165b2f1ab58SBarry Smith 1662205254eSKarl Rupp if (nnz1<nnz2) return -1; 1672205254eSKarl Rupp if (nnz1>nnz2) return 1; 168b2f1ab58SBarry Smith 1694e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 1704e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 1712205254eSKarl Rupp if (icol1[j]< icol2[j]) return -1; 1722205254eSKarl Rupp if (icol1[j]> icol2[j]) return 1; 173b2f1ab58SBarry Smith } 1744e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 1754e1ff37aSBarry Smith for (j=0; j<nnz1; j++) { 1762205254eSKarl Rupp if (icol1[j]-i1< icol2[j]-i2) return -1; 1772205254eSKarl Rupp if (icol1[j]-i1> icol2[j]-i2) return 1; 178b2f1ab58SBarry Smith } 1794e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 1804e1ff37aSBarry Smith for (j=1; j<nnz1; j++) { 1812205254eSKarl Rupp if (icol1[j]-icol1[0] < icol2[j]-icol2[0]) return -1; 1822205254eSKarl Rupp if (icol1[j]-icol1[0] > icol2[j]-icol2[0]) return 1; 183b2f1ab58SBarry Smith } 184b2f1ab58SBarry Smith } 185b2f1ab58SBarry Smith return 0; 186b2f1ab58SBarry Smith } 187b2f1ab58SBarry Smith 188b2f1ab58SBarry Smith /* 189b2f1ab58SBarry Smith spbas_mergesort_icols: 190b2f1ab58SBarry Smith return a sorting of the rows in which identical sparseness patterns are 191b2f1ab58SBarry Smith next to each other 192b2f1ab58SBarry Smith */ 193b2f1ab58SBarry Smith #undef __FUNCT__ 194b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort_icols" 195b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort_icols(PetscInt nrows, PetscInt * irow_in, PetscInt * icol_in,PetscInt col_idx_type, PetscInt *isort) 196b2f1ab58SBarry Smith { 197b2f1ab58SBarry Smith PetscErrorCode ierr; 198c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 199c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 200c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 201c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 202c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 203c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 204c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 205b2f1ab58SBarry Smith 206b2f1ab58SBarry Smith PetscFunctionBegin; 207785e854fSJed Brown ierr = PetscMalloc1(nrows,&ialloc);CHKERRQ(ierr); 208b2f1ab58SBarry Smith 209b2f1ab58SBarry Smith ihlp1 = ialloc; 210b2f1ab58SBarry Smith ihlp2 = isort; 211b2f1ab58SBarry Smith 212c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 2134e1ff37aSBarry Smith for (istep=1; istep<nrows; istep*=2) { 214c328eeadSBarry Smith /* 215c328eeadSBarry Smith Combine sorted parts 216c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 217c328eeadSBarry Smith of ihlp2 and vhlp2 218c328eeadSBarry Smith 219c328eeadSBarry Smith into one sorted part 220c328eeadSBarry Smith istart:istart+2*istep-1 221c328eeadSBarry Smith of ihlp1 and vhlp1 222c328eeadSBarry Smith */ 2234e1ff37aSBarry Smith for (istart=0; istart<nrows; istart+=2*istep) { 224c328eeadSBarry Smith /* Set counters and bound array part endings */ 2252205254eSKarl Rupp i1=istart; i1end = i1+istep; if (i1end>nrows) i1end=nrows; 2262205254eSKarl Rupp i2=istart+istep; i2end = i2+istep; if (i2end>nrows) i2end=nrows; 227b2f1ab58SBarry Smith 228c328eeadSBarry Smith /* Merge the two array parts */ 2294e1ff37aSBarry Smith for (i=istart; i<i2end; i++) { 2304e1ff37aSBarry Smith if (i1<i1end && i2<i2end && spbas_row_order_icol(ihlp2[i1], ihlp2[i2], irow_in, icol_in, col_idx_type) < 0) { 231b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 232b2f1ab58SBarry Smith i1++; 2334e1ff37aSBarry Smith } else if (i2<i2end) { 234b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 235b2f1ab58SBarry Smith i2++; 2364e1ff37aSBarry Smith } else { 237b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 238b2f1ab58SBarry Smith i1++; 239b2f1ab58SBarry Smith } 240b2f1ab58SBarry Smith } 241b2f1ab58SBarry Smith } 242b2f1ab58SBarry Smith 243c328eeadSBarry Smith /* Swap the two array sets */ 244b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 245b2f1ab58SBarry Smith } 246b2f1ab58SBarry Smith 247c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 2484e1ff37aSBarry Smith if (ihlp2 != isort) { 2492205254eSKarl Rupp for (i=0; i<nrows; i++) isort[i] = ihlp2[i]; 250b2f1ab58SBarry Smith } 251b2f1ab58SBarry Smith ierr = PetscFree(ialloc);CHKERRQ(ierr); 252b2f1ab58SBarry Smith PetscFunctionReturn(0); 253b2f1ab58SBarry Smith } 254b2f1ab58SBarry Smith 255b2f1ab58SBarry Smith 256b2f1ab58SBarry Smith 257b2f1ab58SBarry Smith /* 258b2f1ab58SBarry Smith spbas_compress_pattern: 259b2f1ab58SBarry Smith calculate a compressed sparseness pattern for a sparseness pattern 260b2f1ab58SBarry Smith given in compressed row storage. The compressed sparseness pattern may 261b2f1ab58SBarry Smith require (much) less memory. 262b2f1ab58SBarry Smith */ 263b2f1ab58SBarry Smith #undef __FUNCT__ 264b2f1ab58SBarry Smith #define __FUNCT__ "spbas_compress_pattern" 2654e1ff37aSBarry 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) 266b2f1ab58SBarry Smith { 267b2f1ab58SBarry Smith PetscInt nnz = irow_in[nrows]; 268b2f1ab58SBarry Smith long int mem_orig = (nrows + nnz) * sizeof(PetscInt); 269b2f1ab58SBarry Smith long int mem_compressed; 270b2f1ab58SBarry Smith PetscErrorCode ierr; 271b2f1ab58SBarry Smith PetscInt *isort; 272b2f1ab58SBarry Smith PetscInt *icols; 273b2f1ab58SBarry Smith PetscInt row_nnz; 274b2f1ab58SBarry Smith PetscInt *ipoint; 275ace3abfcSBarry Smith PetscBool *used; 276b2f1ab58SBarry Smith PetscInt ptr; 277b2f1ab58SBarry Smith PetscInt i,j; 278ace3abfcSBarry Smith const PetscBool no_values = PETSC_FALSE; 279b2f1ab58SBarry Smith 280b2f1ab58SBarry Smith PetscFunctionBegin; 281c328eeadSBarry Smith /* Allocate the structure of the new matrix */ 282b2f1ab58SBarry Smith B->nrows = nrows; 283b2f1ab58SBarry Smith B->ncols = ncols; 284b2f1ab58SBarry Smith B->nnz = nnz; 285b2f1ab58SBarry Smith B->col_idx_type = col_idx_type; 286b2f1ab58SBarry Smith B->block_data = PETSC_TRUE; 2872205254eSKarl Rupp 288b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(B, no_values);CHKERRQ(ierr); 289b2f1ab58SBarry Smith 290c328eeadSBarry Smith /* When using an offset array, set it */ 2914e1ff37aSBarry Smith if (col_idx_type==SPBAS_OFFSET_ARRAY) { 2922205254eSKarl Rupp for (i=0; i<nrows; i++) B->icol0[i] = icol_in[irow_in[i]]; 293b2f1ab58SBarry Smith } 294b2f1ab58SBarry Smith 295c328eeadSBarry Smith /* Allocate the ordering for the rows */ 296785e854fSJed Brown ierr = PetscMalloc1(nrows,&isort);CHKERRQ(ierr); 297785e854fSJed Brown ierr = PetscMalloc1(nrows,&ipoint);CHKERRQ(ierr); 298785e854fSJed Brown ierr = PetscMalloc1(nrows,&used);CHKERRQ(ierr); 299b2f1ab58SBarry Smith 300c328eeadSBarry Smith /* Initialize the sorting */ 301ace3abfcSBarry Smith ierr = PetscMemzero((void*) used, nrows*sizeof(PetscBool));CHKERRQ(ierr); 3024e1ff37aSBarry Smith for (i = 0; i<nrows; i++) { 303b2f1ab58SBarry Smith B->row_nnz[i] = irow_in[i+1]-irow_in[i]; 304b2f1ab58SBarry Smith isort[i] = i; 305b2f1ab58SBarry Smith ipoint[i] = i; 306b2f1ab58SBarry Smith } 307b2f1ab58SBarry Smith 308c328eeadSBarry Smith /* Sort the rows so that identical columns will be next to each other */ 309b2f1ab58SBarry Smith ierr = spbas_mergesort_icols(nrows, irow_in, icol_in, col_idx_type, isort);CHKERRQ(ierr); 3100298fd71SBarry Smith ierr = PetscInfo(NULL,"Rows have been sorted for patterns\n");CHKERRQ(ierr); 311b2f1ab58SBarry Smith 312c328eeadSBarry Smith /* Replace identical rows with the first one in the list */ 3134e1ff37aSBarry Smith for (i=1; i<nrows; i++) { 3144e1ff37aSBarry Smith if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) { 315b2f1ab58SBarry Smith ipoint[isort[i]] = ipoint[isort[i-1]]; 316b2f1ab58SBarry Smith } 317b2f1ab58SBarry Smith } 318b2f1ab58SBarry Smith 319c328eeadSBarry Smith /* Collect the rows which are used*/ 3202205254eSKarl Rupp for (i=0; i<nrows; i++) used[ipoint[i]] = PETSC_TRUE; 321b2f1ab58SBarry Smith 322c328eeadSBarry Smith /* Calculate needed memory */ 323b2f1ab58SBarry Smith B->n_alloc_icol = 0; 3244e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 3252205254eSKarl Rupp if (used[i]) B->n_alloc_icol += B->row_nnz[i]; 326b2f1ab58SBarry Smith } 327785e854fSJed Brown ierr = PetscMalloc1(B->n_alloc_icol,&B->alloc_icol);CHKERRQ(ierr); 328b2f1ab58SBarry Smith 329c328eeadSBarry Smith /* Fill in the diagonal offsets for the rows which store their own data */ 330b2f1ab58SBarry Smith ptr = 0; 3314e1ff37aSBarry Smith for (i=0; i<B->nrows; i++) { 3324e1ff37aSBarry Smith if (used[i]) { 333b2f1ab58SBarry Smith B->icols[i] = &B->alloc_icol[ptr]; 334b2f1ab58SBarry Smith icols = &icol_in[irow_in[i]]; 335b2f1ab58SBarry Smith row_nnz = B->row_nnz[i]; 3364e1ff37aSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 3374e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 338b2f1ab58SBarry Smith B->icols[i][j] = icols[j]; 339b2f1ab58SBarry Smith } 3404e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 3414e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 342b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-i; 343b2f1ab58SBarry Smith } 3444e1ff37aSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 3454e1ff37aSBarry Smith for (j=0; j<row_nnz; j++) { 346b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-icols[0]; 347b2f1ab58SBarry Smith } 348b2f1ab58SBarry Smith } 349b2f1ab58SBarry Smith ptr += B->row_nnz[i]; 350b2f1ab58SBarry Smith } 351b2f1ab58SBarry Smith } 352b2f1ab58SBarry Smith 353c328eeadSBarry Smith /* Point to the right places for all data */ 3544e1ff37aSBarry Smith for (i=0; i<nrows; i++) { 355b2f1ab58SBarry Smith B->icols[i] = B->icols[ipoint[i]]; 356b2f1ab58SBarry Smith } 3570298fd71SBarry Smith ierr = PetscInfo(NULL,"Row patterns have been compressed\n");CHKERRQ(ierr); 35857622a8eSBarry Smith ierr = PetscInfo1(NULL," (%g nonzeros per row)\n", (double) ((PetscReal) nnz / (PetscReal) nrows));CHKERRQ(ierr); 359b2f1ab58SBarry Smith 360b2f1ab58SBarry Smith ierr=PetscFree(isort);CHKERRQ(ierr); 361b2f1ab58SBarry Smith ierr=PetscFree(used);CHKERRQ(ierr); 362b2f1ab58SBarry Smith ierr=PetscFree(ipoint);CHKERRQ(ierr); 363b2f1ab58SBarry Smith 364b2f1ab58SBarry Smith mem_compressed = spbas_memory_requirement(*B); 3654e1ff37aSBarry Smith *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig; 366b2f1ab58SBarry Smith PetscFunctionReturn(0); 367b2f1ab58SBarry Smith } 368b2f1ab58SBarry Smith 369b2f1ab58SBarry Smith /* 370b2f1ab58SBarry Smith spbas_incomplete_cholesky 371b2f1ab58SBarry Smith Incomplete Cholesky decomposition 372b2f1ab58SBarry Smith */ 373c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas_cholesky.h> 374b2f1ab58SBarry Smith 375b2f1ab58SBarry Smith /* 376b2f1ab58SBarry Smith spbas_delete : de-allocate the arrays owned by this matrix 377b2f1ab58SBarry Smith */ 378b2f1ab58SBarry Smith #undef __FUNCT__ 379b2f1ab58SBarry Smith #define __FUNCT__ "spbas_delete" 380b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix matrix) 381b2f1ab58SBarry Smith { 382b2f1ab58SBarry Smith PetscInt i; 383b2f1ab58SBarry Smith PetscErrorCode ierr; 3849ccc4a9bSBarry Smith 385b2f1ab58SBarry Smith PetscFunctionBegin; 3869ccc4a9bSBarry Smith if (matrix.block_data) { 387cb9801acSJed Brown ierr=PetscFree(matrix.alloc_icol);CHKERRQ(ierr); 388b2f1ab58SBarry Smith if (matrix.values) {ierr=PetscFree(matrix.alloc_val);CHKERRQ(ierr);} 3899ccc4a9bSBarry Smith } else { 3909ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.icols[i]);CHKERRQ(ierr);} 391b2f1ab58SBarry Smith ierr = PetscFree(matrix.icols);CHKERRQ(ierr); 3929ccc4a9bSBarry Smith if (matrix.values) { 3939ccc4a9bSBarry Smith for (i=0; i<matrix.nrows; i++) { ierr=PetscFree(matrix.values[i]);CHKERRQ(ierr);} 394b2f1ab58SBarry Smith } 395b2f1ab58SBarry Smith } 396b2f1ab58SBarry Smith 397b2f1ab58SBarry Smith ierr=PetscFree(matrix.row_nnz);CHKERRQ(ierr); 398b2f1ab58SBarry Smith ierr=PetscFree(matrix.icols);CHKERRQ(ierr); 3999ccc4a9bSBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) {ierr=PetscFree(matrix.icol0);CHKERRQ(ierr);} 400c31cb41cSBarry Smith ierr=PetscFree(matrix.values);CHKERRQ(ierr); 401b2f1ab58SBarry Smith PetscFunctionReturn(0); 402b2f1ab58SBarry Smith } 403b2f1ab58SBarry Smith 404b2f1ab58SBarry Smith /* 405b2f1ab58SBarry Smith spbas_matrix_to_crs: 406b2f1ab58SBarry Smith Convert an spbas_matrix to compessed row storage 407b2f1ab58SBarry Smith */ 408b2f1ab58SBarry Smith #undef __FUNCT__ 409b2f1ab58SBarry Smith #define __FUNCT__ "spbas_matrix_to_crs" 410b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out) 411b2f1ab58SBarry Smith { 412b2f1ab58SBarry Smith PetscInt nrows = matrix_A.nrows; 413b2f1ab58SBarry Smith PetscInt nnz = matrix_A.nnz; 414b2f1ab58SBarry Smith PetscInt i,j,r_nnz,i0; 415b2f1ab58SBarry Smith PetscInt *irow; 416b2f1ab58SBarry Smith PetscInt *icol; 417b2f1ab58SBarry Smith PetscInt *icol_A; 418b2f1ab58SBarry Smith MatScalar *val; 419b2f1ab58SBarry Smith PetscScalar *val_A; 420b2f1ab58SBarry Smith PetscInt col_idx_type = matrix_A.col_idx_type; 421ace3abfcSBarry Smith PetscBool do_values = matrix_A.values ? PETSC_TRUE : PETSC_FALSE; 422b2f1ab58SBarry Smith PetscErrorCode ierr; 423b2f1ab58SBarry Smith 424b2f1ab58SBarry Smith PetscFunctionBegin; 425854ce69bSBarry Smith ierr = PetscMalloc1(nrows+1, &irow);CHKERRQ(ierr); 426854ce69bSBarry Smith ierr = PetscMalloc1(nnz, &icol);CHKERRQ(ierr); 4279ccc4a9bSBarry Smith *icol_out = icol; 4289ccc4a9bSBarry Smith *irow_out = irow; 4299ccc4a9bSBarry Smith if (do_values) { 430854ce69bSBarry Smith ierr = PetscMalloc1(nnz, &val);CHKERRQ(ierr); 431b2f1ab58SBarry Smith *val_out = val; *icol_out = icol; *irow_out=irow; 432b2f1ab58SBarry Smith } 433b2f1ab58SBarry Smith 434b2f1ab58SBarry Smith irow[0]=0; 4359ccc4a9bSBarry Smith for (i=0; i<nrows; i++) { 436b2f1ab58SBarry Smith r_nnz = matrix_A.row_nnz[i]; 437b2f1ab58SBarry Smith i0 = irow[i]; 438b2f1ab58SBarry Smith irow[i+1] = i0 + r_nnz; 439b2f1ab58SBarry Smith icol_A = matrix_A.icols[i]; 440b2f1ab58SBarry Smith 4419ccc4a9bSBarry Smith if (do_values) { 442b2f1ab58SBarry Smith val_A = matrix_A.values[i]; 4439ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 444b2f1ab58SBarry Smith icol[i0+j] = icol_A[j]; 445b2f1ab58SBarry Smith val[i0+j] = val_A[j]; 446b2f1ab58SBarry Smith } 4479ccc4a9bSBarry Smith } else { 4482205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] = icol_A[j]; 449b2f1ab58SBarry Smith } 450b2f1ab58SBarry Smith 4519ccc4a9bSBarry Smith if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 4522205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] += i; 4532205254eSKarl Rupp } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 454b2f1ab58SBarry Smith i0 = matrix_A.icol0[i]; 4552205254eSKarl Rupp for (j=0; j<r_nnz; j++) icol[i0+j] += i0; 456b2f1ab58SBarry Smith } 457b2f1ab58SBarry Smith } 458b2f1ab58SBarry Smith PetscFunctionReturn(0); 459b2f1ab58SBarry Smith } 460b2f1ab58SBarry Smith 461b2f1ab58SBarry Smith 462b2f1ab58SBarry Smith /* 463b2f1ab58SBarry Smith spbas_transpose 464b2f1ab58SBarry Smith return the transpose of a matrix 465b2f1ab58SBarry Smith */ 466b2f1ab58SBarry Smith #undef __FUNCT__ 467b2f1ab58SBarry Smith #define __FUNCT__ "spbas_transpose" 468b2f1ab58SBarry Smith PetscErrorCode spbas_transpose(spbas_matrix in_matrix, spbas_matrix * result) 469b2f1ab58SBarry Smith { 470b2f1ab58SBarry Smith PetscInt col_idx_type = in_matrix.col_idx_type; 471b2f1ab58SBarry Smith PetscInt nnz = in_matrix.nnz; 472b2f1ab58SBarry Smith PetscInt ncols = in_matrix.nrows; 473b2f1ab58SBarry Smith PetscInt nrows = in_matrix.ncols; 474b2f1ab58SBarry Smith PetscInt i,j,k; 475b2f1ab58SBarry Smith PetscInt r_nnz; 476b2f1ab58SBarry Smith PetscInt *irow; 4774efc9174SBarry Smith PetscInt icol0 = 0; 478b2f1ab58SBarry Smith PetscScalar * val; 479b2f1ab58SBarry Smith PetscErrorCode ierr; 4804e1ff37aSBarry Smith 481b2f1ab58SBarry Smith PetscFunctionBegin; 482c328eeadSBarry Smith /* Copy input values */ 483b2f1ab58SBarry Smith result->nrows = nrows; 484b2f1ab58SBarry Smith result->ncols = ncols; 485b2f1ab58SBarry Smith result->nnz = nnz; 486b2f1ab58SBarry Smith result->col_idx_type = SPBAS_COLUMN_NUMBERS; 487b2f1ab58SBarry Smith result->block_data = PETSC_TRUE; 488b2f1ab58SBarry Smith 489c328eeadSBarry Smith /* Allocate sparseness pattern */ 49071d55d18SBarry Smith ierr = spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 491b2f1ab58SBarry Smith 492c328eeadSBarry Smith /* Count the number of nonzeros in each row */ 4932205254eSKarl Rupp for (i = 0; i<nrows; i++) result->row_nnz[i] = 0; 494b2f1ab58SBarry Smith 4959ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 496b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 497b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 4989ccc4a9bSBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) { 4992205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[irow[j]]++; 5009ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { 5012205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[i+irow[j]]++; 5029ccc4a9bSBarry Smith } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { 503b2f1ab58SBarry Smith icol0=in_matrix.icol0[i]; 5042205254eSKarl Rupp for (j=0; j<r_nnz; j++) result->row_nnz[icol0+irow[j]]++; 505b2f1ab58SBarry Smith } 506b2f1ab58SBarry Smith } 507b2f1ab58SBarry Smith 508c328eeadSBarry Smith /* Set the pointers to the data */ 509b2f1ab58SBarry Smith ierr = spbas_allocate_data(result);CHKERRQ(ierr); 510b2f1ab58SBarry Smith 511c328eeadSBarry Smith /* Reset the number of nonzeros in each row */ 5122205254eSKarl Rupp for (i = 0; i<nrows; i++) result->row_nnz[i] = 0; 513b2f1ab58SBarry Smith 514c328eeadSBarry Smith /* Fill the data arrays */ 5159ccc4a9bSBarry Smith if (in_matrix.values) { 5169ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 517b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 518b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 519b2f1ab58SBarry Smith val = in_matrix.values[i]; 520b2f1ab58SBarry Smith 5212205254eSKarl Rupp if (col_idx_type == SPBAS_COLUMN_NUMBERS) icol0 = 0; 5222205254eSKarl Rupp else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0 = i; 5232205254eSKarl Rupp else if (col_idx_type == SPBAS_OFFSET_ARRAY) icol0 = in_matrix.icol0[i]; 5249ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 525b2f1ab58SBarry Smith k = icol0 + irow[j]; 526b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 527b2f1ab58SBarry Smith result->values[k][result->row_nnz[k]] = val[j]; 528b2f1ab58SBarry Smith result->row_nnz[k]++; 529b2f1ab58SBarry Smith } 530b2f1ab58SBarry Smith } 5319ccc4a9bSBarry Smith } else { 5329ccc4a9bSBarry Smith for (i=0; i<ncols; i++) { 533b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 534b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 535b2f1ab58SBarry Smith 5362205254eSKarl Rupp if (col_idx_type == SPBAS_COLUMN_NUMBERS) icol0=0; 5372205254eSKarl Rupp else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) icol0=i; 5382205254eSKarl Rupp else if (col_idx_type == SPBAS_OFFSET_ARRAY) icol0=in_matrix.icol0[i]; 539b2f1ab58SBarry Smith 5409ccc4a9bSBarry Smith for (j=0; j<r_nnz; j++) { 541b2f1ab58SBarry Smith k = icol0 + irow[j]; 542b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 543b2f1ab58SBarry Smith result->row_nnz[k]++; 544b2f1ab58SBarry Smith } 545b2f1ab58SBarry Smith } 546b2f1ab58SBarry Smith } 547b2f1ab58SBarry Smith PetscFunctionReturn(0); 548b2f1ab58SBarry Smith } 549b2f1ab58SBarry Smith 550b2f1ab58SBarry Smith /* 551b2f1ab58SBarry Smith spbas_mergesort 552b2f1ab58SBarry Smith 553b2f1ab58SBarry Smith mergesort for an array of intergers and an array of associated 554b2f1ab58SBarry Smith reals 555b2f1ab58SBarry Smith 556b2f1ab58SBarry Smith on output, icol[0..nnz-1] is increasing; 557b2f1ab58SBarry Smith val[0..nnz-1] has undergone the same permutation as icol 558b2f1ab58SBarry Smith 5590298fd71SBarry Smith NB: val may be NULL: in that case, only the integers are sorted 560b2f1ab58SBarry Smith 561b2f1ab58SBarry Smith */ 562b2f1ab58SBarry Smith #undef __FUNCT__ 563b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort" 564b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val) 565b2f1ab58SBarry Smith { 566c328eeadSBarry Smith PetscInt istep; /* Chunk-sizes of already sorted parts of arrays */ 567c328eeadSBarry Smith PetscInt i, i1, i2; /* Loop counters for (partly) sorted arrays */ 568c328eeadSBarry Smith PetscInt istart, i1end, i2end; /* start of newly sorted array part, end of both parts */ 569c328eeadSBarry Smith PetscInt *ialloc; /* Allocated arrays */ 5700298fd71SBarry Smith PetscScalar *valloc=NULL; 571c328eeadSBarry Smith PetscInt *iswap; /* auxiliary pointers for swapping */ 572b2f1ab58SBarry Smith PetscScalar *vswap; 573c328eeadSBarry Smith PetscInt *ihlp1; /* Pointers to new version of arrays, */ 5740298fd71SBarry Smith PetscScalar *vhlp1=NULL; /* (arrays under construction) */ 575c328eeadSBarry Smith PetscInt *ihlp2; /* Pointers to previous version of arrays, */ 5760298fd71SBarry Smith PetscScalar *vhlp2=NULL; 577b2f1ab58SBarry Smith PetscErrorCode ierr; 578b2f1ab58SBarry Smith 579785e854fSJed Brown ierr = PetscMalloc1(nnz,&ialloc);CHKERRQ(ierr); 580b2f1ab58SBarry Smith ihlp1 = ialloc; 581b2f1ab58SBarry Smith ihlp2 = icol; 582b2f1ab58SBarry Smith 5839ccc4a9bSBarry Smith if (val) { 584785e854fSJed Brown ierr = PetscMalloc1(nnz,&valloc);CHKERRQ(ierr); 585b2f1ab58SBarry Smith vhlp1 = valloc; 586b2f1ab58SBarry Smith vhlp2 = val; 587b2f1ab58SBarry Smith } 588b2f1ab58SBarry Smith 589b2f1ab58SBarry Smith 590c328eeadSBarry Smith /* Sorted array chunks are first 1 long, and increase until they are the complete array */ 5919ccc4a9bSBarry Smith for (istep=1; istep<nnz; istep*=2) { 592c328eeadSBarry Smith /* 593c328eeadSBarry Smith Combine sorted parts 594c328eeadSBarry Smith istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 595c328eeadSBarry Smith of ihlp2 and vhlp2 596c328eeadSBarry Smith 597c328eeadSBarry Smith into one sorted part 598c328eeadSBarry Smith istart:istart+2*istep-1 599c328eeadSBarry Smith of ihlp1 and vhlp1 600c328eeadSBarry Smith */ 6019ccc4a9bSBarry Smith for (istart=0; istart<nnz; istart+=2*istep) { 602c328eeadSBarry Smith /* Set counters and bound array part endings */ 6032205254eSKarl Rupp i1=istart; i1end = i1+istep; if (i1end>nnz) i1end=nnz; 6042205254eSKarl Rupp i2=istart+istep; i2end = i2+istep; if (i2end>nnz) i2end=nnz; 605b2f1ab58SBarry Smith 606c328eeadSBarry Smith /* Merge the two array parts */ 6079ccc4a9bSBarry Smith if (val) { 6089ccc4a9bSBarry Smith for (i=istart; i<i2end; i++) { 6099ccc4a9bSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 610b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 611b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 612b2f1ab58SBarry Smith i1++; 6139ccc4a9bSBarry Smith } else if (i2<i2end) { 614b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 615b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i2]; 616b2f1ab58SBarry Smith i2++; 6179ccc4a9bSBarry Smith } else { 618b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 619b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 620b2f1ab58SBarry Smith i1++; 621b2f1ab58SBarry Smith } 622b2f1ab58SBarry Smith } 6239ccc4a9bSBarry Smith } else { 6246322f4bdSBarry Smith for (i=istart; i<i2end; i++) { 6256322f4bdSBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) { 626b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 627b2f1ab58SBarry Smith i1++; 6286322f4bdSBarry Smith } else if (i2<i2end) { 629b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 630b2f1ab58SBarry Smith i2++; 6316322f4bdSBarry Smith } else { 632b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 633b2f1ab58SBarry Smith i1++; 634b2f1ab58SBarry Smith } 635b2f1ab58SBarry Smith } 636b2f1ab58SBarry Smith } 637b2f1ab58SBarry Smith } 638b2f1ab58SBarry Smith 639c328eeadSBarry Smith /* Swap the two array sets */ 640b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 641b2f1ab58SBarry Smith vswap = vhlp2; vhlp2 = vhlp1; vhlp1 = vswap; 642b2f1ab58SBarry Smith } 643b2f1ab58SBarry Smith 644c328eeadSBarry Smith /* Copy one more time in case the sorted arrays are the temporary ones */ 6456322f4bdSBarry Smith if (ihlp2 != icol) { 6462205254eSKarl Rupp for (i=0; i<nnz; i++) icol[i] = ihlp2[i]; 6476322f4bdSBarry Smith if (val) { 6482205254eSKarl Rupp for (i=0; i<nnz; i++) val[i] = vhlp2[i]; 649b2f1ab58SBarry Smith } 650b2f1ab58SBarry Smith } 651b2f1ab58SBarry Smith 652b2f1ab58SBarry Smith ierr = PetscFree(ialloc);CHKERRQ(ierr); 653b2f1ab58SBarry Smith if (val) {ierr = PetscFree(valloc);CHKERRQ(ierr);} 654b2f1ab58SBarry Smith PetscFunctionReturn(0); 655b2f1ab58SBarry Smith } 656b2f1ab58SBarry Smith 657b2f1ab58SBarry Smith /* 658b2f1ab58SBarry Smith spbas_apply_reordering_rows: 659b2f1ab58SBarry Smith apply the given reordering to the rows: matrix_A = matrix_A(perm,:); 660b2f1ab58SBarry Smith */ 661b2f1ab58SBarry Smith #undef __FUNCT__ 662b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_rows" 663b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_rows(spbas_matrix *matrix_A, const PetscInt *permutation) 664b2f1ab58SBarry Smith { 665b2f1ab58SBarry Smith PetscInt i,j,ip; 666b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 667b2f1ab58SBarry Smith PetscInt * row_nnz; 668b2f1ab58SBarry Smith PetscInt **icols; 669ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 6700298fd71SBarry Smith PetscScalar **vals = NULL; 671b2f1ab58SBarry Smith PetscErrorCode ierr; 672b2f1ab58SBarry Smith 673b2f1ab58SBarry Smith PetscFunctionBegin; 674e32f2f54SBarry 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"); 675b2f1ab58SBarry Smith 6766322f4bdSBarry Smith if (do_values) { 677854ce69bSBarry Smith ierr = PetscMalloc1(nrows, &vals);CHKERRQ(ierr); 678b2f1ab58SBarry Smith } 679854ce69bSBarry Smith ierr = PetscMalloc1(nrows, &row_nnz);CHKERRQ(ierr); 680854ce69bSBarry Smith ierr = PetscMalloc1(nrows, &icols);CHKERRQ(ierr); 681b2f1ab58SBarry Smith 6826322f4bdSBarry Smith for (i=0; i<nrows; i++) { 683b2f1ab58SBarry Smith ip = permutation[i]; 6842205254eSKarl Rupp if (do_values) vals[i] = matrix_A->values[ip]; 685b2f1ab58SBarry Smith icols[i] = matrix_A->icols[ip]; 686b2f1ab58SBarry Smith row_nnz[i] = matrix_A->row_nnz[ip]; 6872205254eSKarl Rupp for (j=0; j<row_nnz[i]; j++) icols[i][j] += ip-i; 688b2f1ab58SBarry Smith } 689b2f1ab58SBarry Smith 690b2f1ab58SBarry Smith if (do_values) { ierr = PetscFree(matrix_A->values);CHKERRQ(ierr);} 691b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->icols);CHKERRQ(ierr); 692b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->row_nnz);CHKERRQ(ierr); 693b2f1ab58SBarry Smith 6942205254eSKarl Rupp if (do_values) matrix_A->values = vals; 695b2f1ab58SBarry Smith matrix_A->icols = icols; 696b2f1ab58SBarry Smith matrix_A->row_nnz = row_nnz; 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; 713ace3abfcSBarry Smith PetscBool do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 7140298fd71SBarry Smith PetscScalar *vals = NULL; 715b2f1ab58SBarry Smith PetscErrorCode ierr; 716b2f1ab58SBarry Smith 717b2f1ab58SBarry Smith PetscFunctionBegin; 718e32f2f54SBarry 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"); 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]; 7232205254eSKarl Rupp 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; 7425fd66863SKarl Rupp 743b2f1ab58SBarry Smith PetscFunctionBegin; 744b2f1ab58SBarry Smith ierr = spbas_apply_reordering_rows(matrix_A, inv_perm);CHKERRQ(ierr); 745b2f1ab58SBarry Smith ierr = spbas_apply_reordering_cols(matrix_A, permutation);CHKERRQ(ierr); 746b2f1ab58SBarry Smith PetscFunctionReturn(0); 747b2f1ab58SBarry Smith } 748b2f1ab58SBarry Smith 749b2f1ab58SBarry Smith #undef __FUNCT__ 750b2f1ab58SBarry Smith #define __FUNCT__ "spbas_pattern_only" 751b2f1ab58SBarry Smith PetscErrorCode spbas_pattern_only(PetscInt nrows, PetscInt ncols, PetscInt *ai, PetscInt *aj, spbas_matrix * result) 752b2f1ab58SBarry Smith { 753b2f1ab58SBarry Smith spbas_matrix retval; 754b2f1ab58SBarry Smith PetscInt i, j, i0, r_nnz; 755b2f1ab58SBarry Smith PetscErrorCode ierr; 756b2f1ab58SBarry Smith 757b2f1ab58SBarry Smith PetscFunctionBegin; 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); 7682205254eSKarl Rupp 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" 791be332245SKarl Rupp PetscErrorCode spbas_mark_row_power(PetscInt *iwork, /* marker-vector */ 792c328eeadSBarry Smith PetscInt row, /* row for which the columns are marked */ 793c328eeadSBarry Smith spbas_matrix * in_matrix, /* matrix for which the power is being calculated */ 794c328eeadSBarry Smith PetscInt marker, /* marker-value: 2^power */ 795c328eeadSBarry Smith PetscInt minmrk, /* lower bound for marked points */ 796c328eeadSBarry Smith PetscInt maxmrk) /* upper bound for marked points */ 797b2f1ab58SBarry Smith { 798b2f1ab58SBarry Smith PetscErrorCode ierr; 799b2f1ab58SBarry Smith PetscInt i,j, nnz; 800b2f1ab58SBarry Smith 801b2f1ab58SBarry Smith PetscFunctionBegin; 802b2f1ab58SBarry Smith nnz = in_matrix->row_nnz[row]; 803b2f1ab58SBarry Smith 804c328eeadSBarry Smith /* For higher powers, call this function recursively */ 8056322f4bdSBarry Smith if (marker>1) { 8066322f4bdSBarry Smith for (i=0; i<nnz; i++) { 807b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 8086322f4bdSBarry Smith if (minmrk<=j && j<maxmrk && iwork[j] < marker) { 80971d55d18SBarry Smith ierr = spbas_mark_row_power(iwork, row + in_matrix->icols[row][i],in_matrix, marker/2,minmrk,maxmrk);CHKERRQ(ierr); 810b2f1ab58SBarry Smith iwork[j] |= marker; 811b2f1ab58SBarry Smith } 812b2f1ab58SBarry Smith } 8136322f4bdSBarry Smith } else { 814c328eeadSBarry Smith /* Mark the columns reached */ 8156322f4bdSBarry Smith for (i=0; i<nnz; i++) { 816b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 8172205254eSKarl Rupp if (minmrk<=j && j<maxmrk) iwork[j] |= 1; 818b2f1ab58SBarry Smith } 819b2f1ab58SBarry Smith } 820b2f1ab58SBarry Smith PetscFunctionReturn(0); 821b2f1ab58SBarry Smith } 822b2f1ab58SBarry Smith 823b2f1ab58SBarry Smith 824b2f1ab58SBarry Smith /* 825b2f1ab58SBarry Smith spbas_power 826b2f1ab58SBarry Smith Calculate sparseness patterns for incomplete Cholesky decompositions 827b2f1ab58SBarry Smith of a given order: (almost) all nonzeros of the matrix^(order+1) which 828b2f1ab58SBarry Smith are inside the band width are found and stored in the output sparseness 829b2f1ab58SBarry Smith pattern. 830b2f1ab58SBarry Smith */ 831b2f1ab58SBarry Smith #undef __FUNCT__ 832b2f1ab58SBarry Smith #define __FUNCT__ "spbas_power" 833b2f1ab58SBarry Smith PetscErrorCode spbas_power(spbas_matrix in_matrix,PetscInt power, spbas_matrix * result) 834b2f1ab58SBarry Smith { 835b2f1ab58SBarry Smith spbas_matrix retval; 836b2f1ab58SBarry Smith PetscInt nrows = in_matrix.nrows; 837b2f1ab58SBarry Smith PetscInt ncols = in_matrix.ncols; 838b2f1ab58SBarry Smith PetscInt i, j, kend; 839b2f1ab58SBarry Smith PetscInt nnz, inz; 840b2f1ab58SBarry Smith PetscInt *iwork; 841b2f1ab58SBarry Smith PetscInt marker; 842b2f1ab58SBarry Smith PetscInt maxmrk=0; 843b2f1ab58SBarry Smith PetscErrorCode ierr; 844b2f1ab58SBarry Smith 845b2f1ab58SBarry Smith PetscFunctionBegin; 846e32f2f54SBarry 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"); 847e32f2f54SBarry Smith if (ncols != nrows) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Dimension error\n"); 848e32f2f54SBarry Smith if (in_matrix.values) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Input array must be sparseness pattern (no values)"); 849e32f2f54SBarry Smith if (power<=0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Power must be 1 or up"); 850b2f1ab58SBarry Smith 851c328eeadSBarry Smith /* Copy input values*/ 852b2f1ab58SBarry Smith retval.nrows = ncols; 853b2f1ab58SBarry Smith retval.ncols = nrows; 854b2f1ab58SBarry Smith retval.nnz = 0; 855b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 856b2f1ab58SBarry Smith retval.block_data = PETSC_FALSE; 857b2f1ab58SBarry Smith 858c328eeadSBarry Smith /* Allocate sparseness pattern */ 85971d55d18SBarry Smith ierr = spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 860b2f1ab58SBarry Smith 861c328eeadSBarry Smith /* Allocate marker array */ 862785e854fSJed Brown ierr = PetscMalloc1(nrows, &iwork);CHKERRQ(ierr); 863b2f1ab58SBarry Smith 864c328eeadSBarry Smith /* Erase the pattern for this row */ 8654e1ff37aSBarry Smith ierr = PetscMemzero((void*) iwork, retval.nrows*sizeof(PetscInt));CHKERRQ(ierr); 866b2f1ab58SBarry Smith 867c328eeadSBarry Smith /* Calculate marker values */ 8682205254eSKarl Rupp marker = 1; for (i=1; i<power; i++) marker*=2; 869b2f1ab58SBarry Smith 8706322f4bdSBarry Smith for (i=0; i<nrows; i++) { 871c328eeadSBarry Smith /* Calculate the pattern for each row */ 872b2f1ab58SBarry Smith 873b2f1ab58SBarry Smith nnz = in_matrix.row_nnz[i]; 874b2f1ab58SBarry Smith kend = i+in_matrix.icols[i][nnz-1]; 8752205254eSKarl Rupp if (maxmrk<=kend) maxmrk=kend+1; 87671d55d18SBarry Smith ierr = spbas_mark_row_power(iwork, i, &in_matrix, marker, i, maxmrk);CHKERRQ(ierr); 877b2f1ab58SBarry Smith 878c328eeadSBarry Smith /* Count the columns*/ 879b2f1ab58SBarry Smith nnz = 0; 8802205254eSKarl Rupp for (j=i; j<maxmrk; j++) nnz+= (iwork[j]!=0); 881b2f1ab58SBarry Smith 882c328eeadSBarry Smith /* Allocate the column indices */ 883b2f1ab58SBarry Smith retval.row_nnz[i] = nnz; 884785e854fSJed Brown ierr = PetscMalloc1(nnz,&retval.icols[i]);CHKERRQ(ierr); 885b2f1ab58SBarry Smith 886c328eeadSBarry Smith /* Administrate the column indices */ 887b2f1ab58SBarry Smith inz = 0; 8886322f4bdSBarry Smith for (j=i; j<maxmrk; j++) { 8896322f4bdSBarry Smith if (iwork[j]) { 890b2f1ab58SBarry Smith retval.icols[i][inz] = j-i; 891b2f1ab58SBarry Smith inz++; 892b2f1ab58SBarry Smith iwork[j]=0; 893b2f1ab58SBarry Smith } 894b2f1ab58SBarry Smith } 895b2f1ab58SBarry Smith retval.nnz += nnz; 896b2f1ab58SBarry Smith }; 897b2f1ab58SBarry Smith ierr = PetscFree(iwork);CHKERRQ(ierr); 898b2f1ab58SBarry Smith *result = retval; 899b2f1ab58SBarry Smith PetscFunctionReturn(0); 900b2f1ab58SBarry Smith } 901b2f1ab58SBarry Smith 902b2f1ab58SBarry Smith 903b2f1ab58SBarry Smith 904b2f1ab58SBarry Smith /* 905b2f1ab58SBarry Smith spbas_keep_upper: 906b2f1ab58SBarry Smith remove the lower part of the matrix: keep the upper part 907b2f1ab58SBarry Smith */ 908b2f1ab58SBarry Smith #undef __FUNCT__ 909b2f1ab58SBarry Smith #define __FUNCT__ "spbas_keep_upper" 910b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper(spbas_matrix * inout_matrix) 911b2f1ab58SBarry Smith { 912b2f1ab58SBarry Smith PetscInt i, j; 913b2f1ab58SBarry Smith PetscInt jstart; 914b2f1ab58SBarry Smith 915b2f1ab58SBarry Smith PetscFunctionBegin; 916e32f2f54SBarry Smith if (inout_matrix->block_data) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Not yet for block data matrices\n"); 9176322f4bdSBarry Smith for (i=0; i<inout_matrix->nrows; i++) { 9186322f4bdSBarry Smith for (jstart=0; (jstart<inout_matrix->row_nnz[i]) && (inout_matrix->icols[i][jstart]<0); jstart++) {} 9196322f4bdSBarry Smith if (jstart>0) { 9206322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9216322f4bdSBarry Smith inout_matrix->icols[i][j] = inout_matrix->icols[i][j+jstart]; 922b2f1ab58SBarry Smith } 923b2f1ab58SBarry Smith 9246322f4bdSBarry Smith if (inout_matrix->values) { 9256322f4bdSBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) { 9266322f4bdSBarry Smith inout_matrix->values[i][j] = inout_matrix->values[i][j+jstart]; 927b2f1ab58SBarry Smith } 928b2f1ab58SBarry Smith } 929b2f1ab58SBarry Smith 930b2f1ab58SBarry Smith inout_matrix->row_nnz[i] -= jstart; 931b2f1ab58SBarry Smith 9326322f4bdSBarry Smith inout_matrix->icols[i] = (PetscInt*) realloc((void*) inout_matrix->icols[i], inout_matrix->row_nnz[i]*sizeof(PetscInt)); 933b2f1ab58SBarry Smith 9346322f4bdSBarry Smith if (inout_matrix->values) { 9356322f4bdSBarry Smith inout_matrix->values[i] = (PetscScalar*) realloc((void*) inout_matrix->values[i], inout_matrix->row_nnz[i]*sizeof(PetscScalar)); 936b2f1ab58SBarry Smith } 937b2f1ab58SBarry Smith inout_matrix->nnz -= jstart; 938b2f1ab58SBarry Smith } 939b2f1ab58SBarry Smith } 940b2f1ab58SBarry Smith PetscFunctionReturn(0); 941b2f1ab58SBarry Smith } 942b2f1ab58SBarry Smith 943