1b2f1ab58SBarry Smith #include "petsc.h" 2b2f1ab58SBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 3b2f1ab58SBarry Smith #include "spbas.h" 4b2f1ab58SBarry Smith 5b2f1ab58SBarry Smith 6b2f1ab58SBarry Smith 7b2f1ab58SBarry Smith 8b2f1ab58SBarry Smith /* 9b2f1ab58SBarry Smith spbas_memory_requirement: 10b2f1ab58SBarry Smith Calculate the number of bytes needed to store tha matrix 11b2f1ab58SBarry Smith */ 12b2f1ab58SBarry Smith #undef __FUNCT__ 13b2f1ab58SBarry Smith #define __FUNCT__ "spbas_memory_requirement" 14b2f1ab58SBarry Smith long int spbas_memory_requirement( spbas_matrix matrix) 15b2f1ab58SBarry Smith { 16b2f1ab58SBarry Smith 17b2f1ab58SBarry Smith 18b2f1ab58SBarry Smith // Requirement for 19b2f1ab58SBarry Smith long int memreq = 6 * sizeof(PetscInt) + // nrows, ncols, nnz, n_alloc_icol, 20b2f1ab58SBarry Smith // n_alloc_val, col_idx_type 21b2f1ab58SBarry Smith sizeof(PetscTruth) + // block_data 22b2f1ab58SBarry Smith sizeof(PetscScalar**) + // values 23b2f1ab58SBarry Smith sizeof(PetscScalar*) + // alloc_val 24b2f1ab58SBarry Smith 2 * sizeof(int**) + // icols, icols0 25b2f1ab58SBarry Smith 2 * sizeof(PetscInt*) + // row_nnz, alloc_icol 26b2f1ab58SBarry Smith matrix.nrows * sizeof(PetscInt) + // row_nnz[*] 27b2f1ab58SBarry Smith matrix.nrows * sizeof(PetscInt*); // icols[*] 28b2f1ab58SBarry Smith 29b2f1ab58SBarry Smith // icol0[*] 30b2f1ab58SBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) 31b2f1ab58SBarry Smith { memreq += matrix.nrows * sizeof(PetscInt); } 32b2f1ab58SBarry Smith 33b2f1ab58SBarry Smith // icols[*][*] 34b2f1ab58SBarry Smith if (matrix.block_data) 35b2f1ab58SBarry Smith { memreq += matrix.n_alloc_icol * sizeof(PetscInt); } 36b2f1ab58SBarry Smith else 37b2f1ab58SBarry Smith { memreq += matrix.nnz * sizeof(PetscInt); } 38b2f1ab58SBarry Smith 39b2f1ab58SBarry Smith if (matrix.values) 40b2f1ab58SBarry Smith { 41b2f1ab58SBarry Smith memreq += matrix.nrows * sizeof(PetscScalar*); // values[*] 42b2f1ab58SBarry Smith // values[*][*] 43b2f1ab58SBarry Smith if (matrix.block_data) 44b2f1ab58SBarry Smith { memreq += matrix.n_alloc_val * sizeof(PetscScalar); } 45b2f1ab58SBarry Smith else 46b2f1ab58SBarry Smith { memreq += matrix.nnz * sizeof(PetscScalar); } 47b2f1ab58SBarry Smith } 48b2f1ab58SBarry Smith 49b2f1ab58SBarry Smith return memreq; 50b2f1ab58SBarry Smith } 51b2f1ab58SBarry Smith 52b2f1ab58SBarry Smith /* 53b2f1ab58SBarry Smith spbas_allocate_pattern: 54b2f1ab58SBarry Smith allocate the pattern arrays row_nnz, icols and optionally values 55b2f1ab58SBarry Smith */ 56b2f1ab58SBarry Smith #undef __FUNCT__ 57b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_pattern" 58b2f1ab58SBarry Smith PetscErrorCode spbas_allocate_pattern( spbas_matrix * result, PetscTruth do_values) 59b2f1ab58SBarry Smith { 60b2f1ab58SBarry Smith PetscErrorCode ierr; 61b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 62b2f1ab58SBarry Smith PetscInt col_idx_type = result->col_idx_type; 63b2f1ab58SBarry Smith 64b2f1ab58SBarry Smith PetscFunctionBegin; 65b2f1ab58SBarry Smith // Allocate sparseness pattern 66b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->row_nnz); CHKERRQ(ierr); 67b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt*),&result->icols); CHKERRQ(ierr); 68b2f1ab58SBarry Smith 69b2f1ab58SBarry Smith // If offsets are given wrt an array, create array 70b2f1ab58SBarry Smith if (col_idx_type == SPBAS_OFFSET_ARRAY) 71b2f1ab58SBarry Smith { 72b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&result->icol0); CHKERRQ(ierr); 73b2f1ab58SBarry Smith } 74b2f1ab58SBarry Smith else 75b2f1ab58SBarry Smith { 76b2f1ab58SBarry Smith result->icol0 = NULL; 77b2f1ab58SBarry Smith } 78b2f1ab58SBarry Smith 79b2f1ab58SBarry Smith // If values are given, allocate values array 80b2f1ab58SBarry Smith if (do_values) 81b2f1ab58SBarry Smith { 82b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscScalar*),&result->values); 83b2f1ab58SBarry Smith CHKERRQ(ierr); 84b2f1ab58SBarry Smith } 85b2f1ab58SBarry Smith else 86b2f1ab58SBarry Smith { 87b2f1ab58SBarry Smith result->values = NULL; 88b2f1ab58SBarry Smith } 89b2f1ab58SBarry Smith 90b2f1ab58SBarry Smith PetscFunctionReturn(0); 91b2f1ab58SBarry Smith } 92b2f1ab58SBarry Smith 93b2f1ab58SBarry Smith /* 94b2f1ab58SBarry Smith spbas_allocate_data: 95b2f1ab58SBarry Smith in case of block_data: 96b2f1ab58SBarry Smith Allocate the data arrays alloc_icol and optionally alloc_val, 97b2f1ab58SBarry Smith set appropriate pointers from icols and values; 98b2f1ab58SBarry Smith in case of !block_data: 99b2f1ab58SBarry Smith Allocate the arrays icols[i] and optionally values[i] 100b2f1ab58SBarry Smith */ 101b2f1ab58SBarry Smith #undef __FUNCT__ 102b2f1ab58SBarry Smith #define __FUNCT__ "spbas_allocate_data" 103b2f1ab58SBarry Smith PetscErrorCode spbas_allocate_data( spbas_matrix * result) 104b2f1ab58SBarry Smith { 105b2f1ab58SBarry Smith PetscInt i; 106b2f1ab58SBarry Smith PetscInt nnz = result->nnz; 107b2f1ab58SBarry Smith PetscInt nrows = result->nrows; 108b2f1ab58SBarry Smith PetscInt r_nnz; 109b2f1ab58SBarry Smith PetscErrorCode ierr; 110b2f1ab58SBarry Smith PetscTruth do_values = (result->values != NULL) ? PETSC_TRUE : PETSC_FALSE; 111b2f1ab58SBarry Smith PetscTruth block_data = result->block_data; 112b2f1ab58SBarry Smith 113b2f1ab58SBarry Smith PetscFunctionBegin; 114b2f1ab58SBarry Smith if (block_data) 115b2f1ab58SBarry Smith { 116b2f1ab58SBarry Smith // Allocate the column number array and point to it 117b2f1ab58SBarry Smith result->n_alloc_icol = nnz; 118b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt), &result->alloc_icol); 119b2f1ab58SBarry Smith CHKERRQ(ierr); 120b2f1ab58SBarry Smith result->icols[0]= result->alloc_icol; 121b2f1ab58SBarry Smith for (i=1; i<nrows; i++) 122b2f1ab58SBarry Smith { 123b2f1ab58SBarry Smith result->icols[i] = result->icols[i-1] + result->row_nnz[i-1]; 124b2f1ab58SBarry Smith } 125b2f1ab58SBarry Smith 126b2f1ab58SBarry Smith // Allocate the value array and point to it 127b2f1ab58SBarry Smith if (do_values) 128b2f1ab58SBarry Smith { 129b2f1ab58SBarry Smith result->n_alloc_val= nnz; 130b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscScalar), &result->alloc_val); 131b2f1ab58SBarry Smith CHKERRQ(ierr); 132b2f1ab58SBarry Smith result->values[0]= result->alloc_val; 133b2f1ab58SBarry Smith for (i=1; i<nrows; i++) 134b2f1ab58SBarry Smith { 135b2f1ab58SBarry Smith result->values[i] = result->values[i-1] + result->row_nnz[i-1]; 136b2f1ab58SBarry Smith } 137b2f1ab58SBarry Smith } 138b2f1ab58SBarry Smith } 139b2f1ab58SBarry Smith else 140b2f1ab58SBarry Smith { 141b2f1ab58SBarry Smith for (i=0; i<nrows; i++) 142b2f1ab58SBarry Smith { 143b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 144b2f1ab58SBarry Smith ierr = PetscMalloc(r_nnz*sizeof(PetscInt), &result->icols[i]); 145b2f1ab58SBarry Smith CHKERRQ(ierr); 146b2f1ab58SBarry Smith } 147b2f1ab58SBarry Smith if (do_values) 148b2f1ab58SBarry Smith { 149b2f1ab58SBarry Smith for (i=0; i<nrows; i++) 150b2f1ab58SBarry Smith { 151b2f1ab58SBarry Smith r_nnz = result->row_nnz[i]; 152b2f1ab58SBarry Smith ierr = PetscMalloc(r_nnz*sizeof(PetscScalar), 153b2f1ab58SBarry Smith &result->values[i]); CHKERRQ(ierr); 154b2f1ab58SBarry Smith } 155b2f1ab58SBarry Smith } 156b2f1ab58SBarry Smith } 157b2f1ab58SBarry Smith 158b2f1ab58SBarry Smith PetscFunctionReturn(0); 159b2f1ab58SBarry Smith } 160b2f1ab58SBarry Smith 161b2f1ab58SBarry Smith /* 162b2f1ab58SBarry Smith spbas_row_order_icol 163b2f1ab58SBarry Smith determine if row i1 should come 164b2f1ab58SBarry Smith + before row i2 in the sorted rows (return -1), 165b2f1ab58SBarry Smith + after (return 1) 166b2f1ab58SBarry Smith + is identical (return 0). 167b2f1ab58SBarry Smith */ 168b2f1ab58SBarry Smith #undef __FUNCT__ 169b2f1ab58SBarry Smith #define __FUNCT__ "spbas_row_order_icol" 170b2f1ab58SBarry Smith int spbas_row_order_icol(PetscInt i1, PetscInt i2, PetscInt *irow_in, PetscInt *icol_in,PetscInt col_idx_type) 171b2f1ab58SBarry Smith { 172b2f1ab58SBarry Smith PetscInt j; 173b2f1ab58SBarry Smith PetscInt nnz1 = irow_in[i1+1] - irow_in[i1]; 174b2f1ab58SBarry Smith PetscInt nnz2 = irow_in[i2+1] - irow_in[i2]; 175b2f1ab58SBarry Smith 176b2f1ab58SBarry Smith PetscInt * icol1 = &icol_in[irow_in[i1]]; 177b2f1ab58SBarry Smith PetscInt * icol2 = &icol_in[irow_in[i2]]; 178b2f1ab58SBarry Smith 179b2f1ab58SBarry Smith if (nnz1<nnz2) {return -1;} 180b2f1ab58SBarry Smith if (nnz1>nnz2) {return 1;} 181b2f1ab58SBarry Smith 182b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) 183b2f1ab58SBarry Smith { 184b2f1ab58SBarry Smith for (j=0; j<nnz1; j++) 185b2f1ab58SBarry Smith { 186b2f1ab58SBarry Smith if (icol1[j]< icol2[j]) {return -1;} 187b2f1ab58SBarry Smith if (icol1[j]> icol2[j]) {return 1;} 188b2f1ab58SBarry Smith } 189b2f1ab58SBarry Smith } 190b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) 191b2f1ab58SBarry Smith { 192b2f1ab58SBarry Smith for (j=0; j<nnz1; j++) 193b2f1ab58SBarry Smith { 194b2f1ab58SBarry Smith if (icol1[j]-i1< icol2[j]-i2) {return -1;} 195b2f1ab58SBarry Smith if (icol1[j]-i1> icol2[j]-i2) {return 1;} 196b2f1ab58SBarry Smith } 197b2f1ab58SBarry Smith } 198b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) 199b2f1ab58SBarry Smith { 200b2f1ab58SBarry Smith for (j=1; j<nnz1; j++) 201b2f1ab58SBarry Smith { 202b2f1ab58SBarry Smith if (icol1[j]-icol1[0] < icol2[j]-icol2[0]) {return -1;} 203b2f1ab58SBarry Smith if (icol1[j]-icol1[0] > icol2[j]-icol2[0]) {return 1;} 204b2f1ab58SBarry Smith } 205b2f1ab58SBarry Smith } 206b2f1ab58SBarry Smith return 0; 207b2f1ab58SBarry Smith } 208b2f1ab58SBarry Smith 209b2f1ab58SBarry Smith /* 210b2f1ab58SBarry Smith spbas_mergesort_icols: 211b2f1ab58SBarry Smith return a sorting of the rows in which identical sparseness patterns are 212b2f1ab58SBarry Smith next to each other 213b2f1ab58SBarry Smith */ 214b2f1ab58SBarry Smith #undef __FUNCT__ 215b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort_icols" 216b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort_icols( PetscInt nrows, PetscInt * irow_in, PetscInt * icol_in,PetscInt col_idx_type, PetscInt *isort) 217b2f1ab58SBarry Smith { 218b2f1ab58SBarry Smith PetscErrorCode ierr; 219b2f1ab58SBarry Smith PetscInt istep; // Chunk-sizes of already sorted parts of arrays 220b2f1ab58SBarry Smith 221b2f1ab58SBarry Smith PetscInt i, i1, i2; // Loop counters for (partly) sorted arrays 222b2f1ab58SBarry Smith 223b2f1ab58SBarry Smith PetscInt istart, i1end, i2end; // start of newly sorted array part, end of both 224b2f1ab58SBarry Smith // parts 225b2f1ab58SBarry Smith 226b2f1ab58SBarry Smith PetscInt *ialloc; // Allocated arrays 227b2f1ab58SBarry Smith PetscInt *iswap; // auxiliary pointers for swapping 228b2f1ab58SBarry Smith PetscInt *ihlp1; // Pointers to new version of arrays, 229b2f1ab58SBarry Smith PetscInt *ihlp2; // Pointers to previous version of arrays, 230b2f1ab58SBarry Smith 231b2f1ab58SBarry Smith PetscFunctionBegin; 232b2f1ab58SBarry Smith 233b2f1ab58SBarry Smith // Create arrays 234b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&ialloc); CHKERRQ(ierr); 235b2f1ab58SBarry Smith 236b2f1ab58SBarry Smith ihlp1 = ialloc; 237b2f1ab58SBarry Smith ihlp2 = isort; 238b2f1ab58SBarry Smith 239b2f1ab58SBarry Smith // Sorted array chunks are first 1 long, and increase until they are the 240b2f1ab58SBarry Smith // complete array 241b2f1ab58SBarry Smith for (istep=1; istep<nrows; istep*=2) 242b2f1ab58SBarry Smith { 243b2f1ab58SBarry Smith // 244b2f1ab58SBarry Smith // Combine sorted parts 245b2f1ab58SBarry Smith // istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 246b2f1ab58SBarry Smith // of ihlp2 and vhlp2 247b2f1ab58SBarry Smith // 248b2f1ab58SBarry Smith // into one sorted part 249b2f1ab58SBarry Smith // istart:istart+2*istep-1 250b2f1ab58SBarry Smith // of ihlp1 and vhlp1 251b2f1ab58SBarry Smith // 252b2f1ab58SBarry Smith for (istart=0; istart<nrows; istart+=2*istep) 253b2f1ab58SBarry Smith { 254b2f1ab58SBarry Smith 255b2f1ab58SBarry Smith // Set counters and bound array part endings 256b2f1ab58SBarry Smith i1=istart; i1end = i1+istep; if (i1end>nrows) {i1end=nrows;} 257b2f1ab58SBarry Smith i2=istart+istep; i2end = i2+istep; if (i2end>nrows) {i2end=nrows;} 258b2f1ab58SBarry Smith 259b2f1ab58SBarry Smith // Merge the two array parts 260b2f1ab58SBarry Smith for (i=istart; i<i2end; i++) 261b2f1ab58SBarry Smith { 262b2f1ab58SBarry Smith if (i1<i1end && i2<i2end && 263b2f1ab58SBarry Smith spbas_row_order_icol( ihlp2[i1], ihlp2[i2], 264b2f1ab58SBarry Smith irow_in, icol_in, col_idx_type) < 0) 265b2f1ab58SBarry Smith { 266b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 267b2f1ab58SBarry Smith i1++; 268b2f1ab58SBarry Smith } 269b2f1ab58SBarry Smith else if (i2<i2end ) 270b2f1ab58SBarry Smith { 271b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 272b2f1ab58SBarry Smith i2++; 273b2f1ab58SBarry Smith } 274b2f1ab58SBarry Smith else 275b2f1ab58SBarry Smith { 276b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 277b2f1ab58SBarry Smith i1++; 278b2f1ab58SBarry Smith } 279b2f1ab58SBarry Smith } 280b2f1ab58SBarry Smith } 281b2f1ab58SBarry Smith 282b2f1ab58SBarry Smith // Swap the two array sets 283b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 284b2f1ab58SBarry Smith } 285b2f1ab58SBarry Smith 286b2f1ab58SBarry Smith // Copy one more time in case the sorted arrays are the temporary ones 287b2f1ab58SBarry Smith if (ihlp2 != isort) 288b2f1ab58SBarry Smith { 289b2f1ab58SBarry Smith for (i=0; i<nrows; i++) { isort[i] = ihlp2[i]; } 290b2f1ab58SBarry Smith } 291b2f1ab58SBarry Smith 292b2f1ab58SBarry Smith // Remove help arrays 293b2f1ab58SBarry Smith ierr = PetscFree(ialloc); CHKERRQ(ierr); 294b2f1ab58SBarry Smith PetscFunctionReturn(0); 295b2f1ab58SBarry Smith } 296b2f1ab58SBarry Smith 297b2f1ab58SBarry Smith 298b2f1ab58SBarry Smith 299b2f1ab58SBarry Smith /* 300b2f1ab58SBarry Smith spbas_compress_pattern: 301b2f1ab58SBarry Smith calculate a compressed sparseness pattern for a sparseness pattern 302b2f1ab58SBarry Smith given in compressed row storage. The compressed sparseness pattern may 303b2f1ab58SBarry Smith require (much) less memory. 304b2f1ab58SBarry Smith */ 305b2f1ab58SBarry Smith #undef __FUNCT__ 306b2f1ab58SBarry Smith #define __FUNCT__ "spbas_compress_pattern" 307b2f1ab58SBarry Smith PetscErrorCode spbas_compress_pattern(PetscInt *irow_in, PetscInt *icol_in, PetscInt nrows, PetscInt ncols, 308*9767453cSBarry Smith PetscInt col_idx_type, spbas_matrix *B,PetscReal *mem_reduction) 309b2f1ab58SBarry Smith { 310b2f1ab58SBarry Smith PetscInt nnz = irow_in[nrows]; 311b2f1ab58SBarry Smith long int mem_orig = (nrows + nnz) * sizeof(PetscInt); 312b2f1ab58SBarry Smith long int mem_compressed; 313b2f1ab58SBarry Smith PetscErrorCode ierr; 314b2f1ab58SBarry Smith PetscInt *isort; 315b2f1ab58SBarry Smith PetscInt *icols; 316b2f1ab58SBarry Smith PetscInt row_nnz; 317b2f1ab58SBarry Smith PetscInt *ipoint; 318b2f1ab58SBarry Smith PetscTruth *used; 319b2f1ab58SBarry Smith PetscInt ptr; 320b2f1ab58SBarry Smith PetscInt i,j; 321b2f1ab58SBarry Smith 322b2f1ab58SBarry Smith const PetscTruth no_values = PETSC_FALSE; 323b2f1ab58SBarry Smith 324b2f1ab58SBarry Smith PetscFunctionBegin; 325b2f1ab58SBarry Smith 326b2f1ab58SBarry Smith // Allocate the structure of the new matrix 327b2f1ab58SBarry Smith B->nrows = nrows; 328b2f1ab58SBarry Smith B->ncols = ncols; 329b2f1ab58SBarry Smith B->nnz = nnz; 330b2f1ab58SBarry Smith B->col_idx_type= col_idx_type; 331b2f1ab58SBarry Smith B->block_data = PETSC_TRUE; 332b2f1ab58SBarry Smith ierr = spbas_allocate_pattern( B, no_values); CHKERRQ(ierr); 333b2f1ab58SBarry Smith 334b2f1ab58SBarry Smith // When using an offset array, set it 335b2f1ab58SBarry Smith if (col_idx_type==SPBAS_OFFSET_ARRAY) 336b2f1ab58SBarry Smith { 337b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {B->icol0[i] = icol_in[irow_in[i]];} 338b2f1ab58SBarry Smith } 339b2f1ab58SBarry Smith 340b2f1ab58SBarry Smith // Allocate the ordering for the rows 341b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&isort); CHKERRQ(ierr); 342b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscInt),&ipoint); CHKERRQ(ierr); 343b2f1ab58SBarry Smith ierr = PetscMalloc(nrows*sizeof(PetscTruth),&used); CHKERRQ(ierr); 344b2f1ab58SBarry Smith 345b2f1ab58SBarry Smith // Initialize the sorting 346b2f1ab58SBarry Smith memset((void*) used, 0, nrows*sizeof(PetscTruth)); 347b2f1ab58SBarry Smith for (i = 0; i<nrows; i++) 348b2f1ab58SBarry Smith { 349b2f1ab58SBarry Smith B->row_nnz[i] = irow_in[i+1]-irow_in[i]; 350b2f1ab58SBarry Smith isort[i] = i; 351b2f1ab58SBarry Smith ipoint[i]= i; 352b2f1ab58SBarry Smith } 353b2f1ab58SBarry Smith 354b2f1ab58SBarry Smith // Sort the rows so that identical columns will be next to each other 355b2f1ab58SBarry Smith ierr = spbas_mergesort_icols( nrows, irow_in, icol_in, col_idx_type, isort); CHKERRQ(ierr); 356b2f1ab58SBarry Smith printf("Rows have been sorted for patterns\n"); 357b2f1ab58SBarry Smith 358b2f1ab58SBarry Smith // Replace identical rows with the first one in the list 359b2f1ab58SBarry Smith for (i=1; i<nrows; i++) 360b2f1ab58SBarry Smith { 361b2f1ab58SBarry Smith if (spbas_row_order_icol(isort[i-1], isort[i], 362b2f1ab58SBarry Smith irow_in, icol_in, col_idx_type) == 0) 363b2f1ab58SBarry Smith { 364b2f1ab58SBarry Smith ipoint[isort[i]] = ipoint[isort[i-1]]; 365b2f1ab58SBarry Smith } 366b2f1ab58SBarry Smith } 367b2f1ab58SBarry Smith 368b2f1ab58SBarry Smith // Collect the rows which are used 369b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {used[ipoint[i]] = PETSC_TRUE;} 370b2f1ab58SBarry Smith 371b2f1ab58SBarry Smith // Calculate needed memory 372b2f1ab58SBarry Smith B->n_alloc_icol = 0; 373b2f1ab58SBarry Smith for (i=0; i<nrows; i++) 374b2f1ab58SBarry Smith { 375b2f1ab58SBarry Smith if (used[i]) {B->n_alloc_icol += B->row_nnz[i];} 376b2f1ab58SBarry Smith } 377b2f1ab58SBarry Smith 378b2f1ab58SBarry Smith ierr = PetscMalloc(B->n_alloc_icol*sizeof(PetscInt),&B->alloc_icol); 379b2f1ab58SBarry Smith CHKERRQ(ierr); 380b2f1ab58SBarry Smith 381b2f1ab58SBarry Smith // Fill in the diagonal offsets for the rows which store their own data 382b2f1ab58SBarry Smith ptr = 0; 383b2f1ab58SBarry Smith for (i=0; i<B->nrows; i++) 384b2f1ab58SBarry Smith { 385b2f1ab58SBarry Smith if (used[i]) 386b2f1ab58SBarry Smith { 387b2f1ab58SBarry Smith B->icols[i] = &B->alloc_icol[ptr]; 388b2f1ab58SBarry Smith icols = &icol_in[irow_in[i]]; 389b2f1ab58SBarry Smith row_nnz = B->row_nnz[i]; 390b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) 391b2f1ab58SBarry Smith { 392b2f1ab58SBarry Smith for (j=0; j<row_nnz; j++) 393b2f1ab58SBarry Smith { 394b2f1ab58SBarry Smith B->icols[i][j] = icols[j]; 395b2f1ab58SBarry Smith } 396b2f1ab58SBarry Smith } 397b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) 398b2f1ab58SBarry Smith { 399b2f1ab58SBarry Smith for (j=0; j<row_nnz; j++) 400b2f1ab58SBarry Smith { 401b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-i; 402b2f1ab58SBarry Smith } 403b2f1ab58SBarry Smith } 404b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) 405b2f1ab58SBarry Smith { 406b2f1ab58SBarry Smith for (j=0; j<row_nnz; j++) 407b2f1ab58SBarry Smith { 408b2f1ab58SBarry Smith B->icols[i][j] = icols[j]-icols[0]; 409b2f1ab58SBarry Smith } 410b2f1ab58SBarry Smith } 411b2f1ab58SBarry Smith ptr += B->row_nnz[i]; 412b2f1ab58SBarry Smith } 413b2f1ab58SBarry Smith } 414b2f1ab58SBarry Smith 415b2f1ab58SBarry Smith // Point to the right places for all data 416b2f1ab58SBarry Smith for (i=0; i<nrows; i++) 417b2f1ab58SBarry Smith { 418b2f1ab58SBarry Smith B->icols[i] = B->icols[ipoint[i]]; 419b2f1ab58SBarry Smith } 420b2f1ab58SBarry Smith printf("Row patterns have been compressed\n"); 421d36aa34eSBarry Smith printf(" (%6.2f nonzeros per row)\n", (PetscReal) nnz / (PetscReal) nrows); 422b2f1ab58SBarry Smith 423b2f1ab58SBarry Smith 424b2f1ab58SBarry Smith // Clean up 425b2f1ab58SBarry Smith ierr=PetscFree(isort); CHKERRQ(ierr); 426b2f1ab58SBarry Smith ierr=PetscFree(used); CHKERRQ(ierr); 427b2f1ab58SBarry Smith ierr=PetscFree(ipoint); CHKERRQ(ierr); 428b2f1ab58SBarry Smith 429b2f1ab58SBarry Smith mem_compressed = spbas_memory_requirement( *B ); 430*9767453cSBarry Smith *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ 431*9767453cSBarry Smith (PetscReal) mem_orig; 432b2f1ab58SBarry Smith PetscFunctionReturn(0); 433b2f1ab58SBarry Smith } 434b2f1ab58SBarry Smith 435b2f1ab58SBarry Smith /* 436b2f1ab58SBarry Smith spbas_incomplete_cholesky 437b2f1ab58SBarry Smith Incomplete Cholesky decomposition 438b2f1ab58SBarry Smith */ 439b2f1ab58SBarry Smith #include "spbas_cholesky.h" 440b2f1ab58SBarry Smith 441b2f1ab58SBarry Smith /* 442b2f1ab58SBarry Smith spbas_delete : de-allocate the arrays owned by this matrix 443b2f1ab58SBarry Smith */ 444b2f1ab58SBarry Smith #undef __FUNCT__ 445b2f1ab58SBarry Smith #define __FUNCT__ "spbas_delete" 446b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix matrix) 447b2f1ab58SBarry Smith { 448b2f1ab58SBarry Smith PetscInt i; 449b2f1ab58SBarry Smith PetscErrorCode ierr; 450b2f1ab58SBarry Smith PetscFunctionBegin; 451b2f1ab58SBarry Smith if (matrix.block_data) 452b2f1ab58SBarry Smith { 453b2f1ab58SBarry Smith ierr=PetscFree(matrix.alloc_icol); CHKERRQ(ierr) 454b2f1ab58SBarry Smith if (matrix.values){ierr=PetscFree(matrix.alloc_val);CHKERRQ(ierr);} 455b2f1ab58SBarry Smith } 456b2f1ab58SBarry Smith else 457b2f1ab58SBarry Smith { 458b2f1ab58SBarry Smith for (i=0; i<matrix.nrows; i++) 459b2f1ab58SBarry Smith { ierr=PetscFree(matrix.icols[i]); CHKERRQ(ierr);} 460b2f1ab58SBarry Smith ierr = PetscFree(matrix.icols);CHKERRQ(ierr); 461b2f1ab58SBarry Smith if (matrix.values) 462b2f1ab58SBarry Smith { 463b2f1ab58SBarry Smith for (i=0; i<matrix.nrows; i++) 464b2f1ab58SBarry Smith { ierr=PetscFree(matrix.values[i]); CHKERRQ(ierr);} 465b2f1ab58SBarry Smith } 466b2f1ab58SBarry Smith } 467b2f1ab58SBarry Smith 468b2f1ab58SBarry Smith ierr=PetscFree(matrix.row_nnz); CHKERRQ(ierr); 469b2f1ab58SBarry Smith ierr=PetscFree(matrix.icols); CHKERRQ(ierr); 470b2f1ab58SBarry Smith if (matrix.col_idx_type == SPBAS_OFFSET_ARRAY) 471b2f1ab58SBarry Smith {ierr=PetscFree(matrix.icol0);CHKERRQ(ierr);} 472b2f1ab58SBarry Smith if (matrix.values) 473b2f1ab58SBarry Smith {ierr=PetscFree(matrix.values);CHKERRQ(ierr);} 474b2f1ab58SBarry Smith 475b2f1ab58SBarry Smith PetscFunctionReturn(0); 476b2f1ab58SBarry Smith } 477b2f1ab58SBarry Smith 478b2f1ab58SBarry Smith /* 479b2f1ab58SBarry Smith spbas_matrix_to_crs: 480b2f1ab58SBarry Smith Convert an spbas_matrix to compessed row storage 481b2f1ab58SBarry Smith */ 482b2f1ab58SBarry Smith #undef __FUNCT__ 483b2f1ab58SBarry Smith #define __FUNCT__ "spbas_matrix_to_crs" 484b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix matrix_A,MatScalar **val_out, PetscInt **irow_out, PetscInt **icol_out) 485b2f1ab58SBarry Smith { 486b2f1ab58SBarry Smith PetscInt nrows = matrix_A.nrows; 487b2f1ab58SBarry Smith PetscInt nnz = matrix_A.nnz; 488b2f1ab58SBarry Smith PetscInt i,j,r_nnz,i0; 489b2f1ab58SBarry Smith PetscInt *irow; 490b2f1ab58SBarry Smith PetscInt *icol; 491b2f1ab58SBarry Smith PetscInt *icol_A; 492b2f1ab58SBarry Smith MatScalar *val; 493b2f1ab58SBarry Smith PetscScalar *val_A; 494b2f1ab58SBarry Smith PetscInt col_idx_type = matrix_A.col_idx_type; 495d36aa34eSBarry Smith PetscTruth do_values = matrix_A.values ? PETSC_TRUE : PETSC_FALSE; 496b2f1ab58SBarry Smith PetscErrorCode ierr; 497b2f1ab58SBarry Smith 498b2f1ab58SBarry Smith PetscFunctionBegin; 499b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt) * (nrows+1), &irow); CHKERRQ(ierr); 500b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt) * nnz, &icol); CHKERRQ(ierr); 501b2f1ab58SBarry Smith *icol_out = icol; *irow_out=irow; 502b2f1ab58SBarry Smith if (do_values) 503b2f1ab58SBarry Smith { 504b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(MatScalar) * nnz, &val); CHKERRQ(ierr); 505b2f1ab58SBarry Smith *val_out = val; *icol_out = icol; *irow_out=irow; 506b2f1ab58SBarry Smith } 507b2f1ab58SBarry Smith 508b2f1ab58SBarry Smith irow[0]=0; 509b2f1ab58SBarry Smith for (i=0; i<nrows; i++) 510b2f1ab58SBarry Smith { 511b2f1ab58SBarry Smith r_nnz = matrix_A.row_nnz[i]; 512b2f1ab58SBarry Smith i0 = irow[i]; 513b2f1ab58SBarry Smith irow[i+1] = i0 + r_nnz; 514b2f1ab58SBarry Smith icol_A = matrix_A.icols[i]; 515b2f1ab58SBarry Smith 516b2f1ab58SBarry Smith if (do_values) 517b2f1ab58SBarry Smith { 518b2f1ab58SBarry Smith val_A = matrix_A.values[i]; 519b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) 520b2f1ab58SBarry Smith { 521b2f1ab58SBarry Smith icol[i0+j] = icol_A[j]; 522b2f1ab58SBarry Smith val[i0+j] = val_A[j]; 523b2f1ab58SBarry Smith } 524b2f1ab58SBarry Smith } 525b2f1ab58SBarry Smith else 526b2f1ab58SBarry Smith { 527b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] = icol_A[j]; } 528b2f1ab58SBarry Smith } 529b2f1ab58SBarry Smith 530b2f1ab58SBarry Smith if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) 531b2f1ab58SBarry Smith { 532b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] += i; } 533b2f1ab58SBarry Smith } 534b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) 535b2f1ab58SBarry Smith { 536b2f1ab58SBarry Smith i0 = matrix_A.icol0[i]; 537b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { icol[i0+j] += i0; } 538b2f1ab58SBarry Smith } 539b2f1ab58SBarry Smith 540b2f1ab58SBarry Smith } 541b2f1ab58SBarry Smith 542b2f1ab58SBarry Smith PetscFunctionReturn(0); 543b2f1ab58SBarry Smith } 544b2f1ab58SBarry Smith 545b2f1ab58SBarry Smith 546b2f1ab58SBarry Smith /* 547b2f1ab58SBarry Smith spbas_transpose 548b2f1ab58SBarry Smith return the transpose of a matrix 549b2f1ab58SBarry Smith */ 550b2f1ab58SBarry Smith #undef __FUNCT__ 551b2f1ab58SBarry Smith #define __FUNCT__ "spbas_transpose" 552b2f1ab58SBarry Smith PetscErrorCode spbas_transpose( spbas_matrix in_matrix, spbas_matrix * result) 553b2f1ab58SBarry Smith { 554b2f1ab58SBarry Smith PetscInt col_idx_type = in_matrix.col_idx_type; 555b2f1ab58SBarry Smith PetscInt nnz = in_matrix.nnz; 556b2f1ab58SBarry Smith PetscInt ncols = in_matrix.nrows; 557b2f1ab58SBarry Smith PetscInt nrows = in_matrix.ncols; 558b2f1ab58SBarry Smith PetscInt i,j,k; 559b2f1ab58SBarry Smith PetscInt r_nnz; 560b2f1ab58SBarry Smith PetscInt *irow; 561b2f1ab58SBarry Smith PetscInt icol0; 562b2f1ab58SBarry Smith PetscScalar * val; 563b2f1ab58SBarry Smith PetscErrorCode ierr; 564b2f1ab58SBarry Smith 565b2f1ab58SBarry Smith PetscFunctionBegin; 566b2f1ab58SBarry Smith 567b2f1ab58SBarry Smith // Copy input values 568b2f1ab58SBarry Smith result->nrows = nrows; 569b2f1ab58SBarry Smith result->ncols = ncols; 570b2f1ab58SBarry Smith result->nnz = nnz; 571b2f1ab58SBarry Smith result->col_idx_type = SPBAS_COLUMN_NUMBERS; 572b2f1ab58SBarry Smith result->block_data = PETSC_TRUE; 573b2f1ab58SBarry Smith 574b2f1ab58SBarry Smith // Allocate sparseness pattern 575d36aa34eSBarry Smith ierr = spbas_allocate_pattern(result, in_matrix.values ? PETSC_TRUE : PETSC_FALSE); 576b2f1ab58SBarry Smith CHKERRQ(ierr); 577b2f1ab58SBarry Smith 578b2f1ab58SBarry Smith // Count the number of nonzeros in each row 579b2f1ab58SBarry Smith for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; } 580b2f1ab58SBarry Smith 581b2f1ab58SBarry Smith for (i=0; i<ncols; i++) 582b2f1ab58SBarry Smith { 583b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 584b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 585b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) 586b2f1ab58SBarry Smith { 587b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[irow[j]]++; } 588b2f1ab58SBarry Smith } 589b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) 590b2f1ab58SBarry Smith { 591b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[i+irow[j]]++; } 592b2f1ab58SBarry Smith } 593b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) 594b2f1ab58SBarry Smith { 595b2f1ab58SBarry Smith icol0=in_matrix.icol0[i]; 596b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) { result->row_nnz[icol0+irow[j]]++; } 597b2f1ab58SBarry Smith } 598b2f1ab58SBarry Smith } 599b2f1ab58SBarry Smith 600b2f1ab58SBarry Smith // Set the pointers to the data 601b2f1ab58SBarry Smith ierr = spbas_allocate_data(result); CHKERRQ(ierr); 602b2f1ab58SBarry Smith 603b2f1ab58SBarry Smith // Reset the number of nonzeros in each row 604b2f1ab58SBarry Smith for (i = 0; i<nrows; i++) { result->row_nnz[i] = 0; } 605b2f1ab58SBarry Smith 606b2f1ab58SBarry Smith // Fill the data arrays 607b2f1ab58SBarry Smith if (in_matrix.values) 608b2f1ab58SBarry Smith { 609b2f1ab58SBarry Smith for (i=0; i<ncols; i++) 610b2f1ab58SBarry Smith { 611b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 612b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 613b2f1ab58SBarry Smith val = in_matrix.values[i]; 614b2f1ab58SBarry Smith 615b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) {icol0=0;} 616b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;} 617b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) 618b2f1ab58SBarry Smith {icol0=in_matrix.icol0[i];} 619b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) 620b2f1ab58SBarry Smith { 621b2f1ab58SBarry Smith k = icol0 + irow[j]; 622b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 623b2f1ab58SBarry Smith result->values[k][result->row_nnz[k]] = val[j]; 624b2f1ab58SBarry Smith result->row_nnz[k]++; 625b2f1ab58SBarry Smith } 626b2f1ab58SBarry Smith } 627b2f1ab58SBarry Smith } 628b2f1ab58SBarry Smith else 629b2f1ab58SBarry Smith { 630b2f1ab58SBarry Smith for (i=0; i<ncols; i++) 631b2f1ab58SBarry Smith { 632b2f1ab58SBarry Smith r_nnz = in_matrix.row_nnz[i]; 633b2f1ab58SBarry Smith irow = in_matrix.icols[i]; 634b2f1ab58SBarry Smith 635b2f1ab58SBarry Smith if (col_idx_type == SPBAS_COLUMN_NUMBERS) {icol0=0;} 636b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) {icol0=i;} 637b2f1ab58SBarry Smith else if (col_idx_type == SPBAS_OFFSET_ARRAY) 638b2f1ab58SBarry Smith {icol0=in_matrix.icol0[i];} 639b2f1ab58SBarry Smith 640b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) 641b2f1ab58SBarry Smith { 642b2f1ab58SBarry Smith k = icol0 + irow[j]; 643b2f1ab58SBarry Smith result->icols[k][result->row_nnz[k]] = i; 644b2f1ab58SBarry Smith result->row_nnz[k]++; 645b2f1ab58SBarry Smith } 646b2f1ab58SBarry Smith } 647b2f1ab58SBarry Smith } 648b2f1ab58SBarry Smith 649b2f1ab58SBarry Smith // Set return value 650b2f1ab58SBarry Smith PetscFunctionReturn(0); 651b2f1ab58SBarry Smith } 652b2f1ab58SBarry Smith 653b2f1ab58SBarry Smith /* 654b2f1ab58SBarry Smith spbas_mergesort 655b2f1ab58SBarry Smith 656b2f1ab58SBarry Smith mergesort for an array of intergers and an array of associated 657b2f1ab58SBarry Smith reals 658b2f1ab58SBarry Smith 659b2f1ab58SBarry Smith on output, icol[0..nnz-1] is increasing; 660b2f1ab58SBarry Smith val[0..nnz-1] has undergone the same permutation as icol 661b2f1ab58SBarry Smith 662b2f1ab58SBarry Smith NB: val may be NULL: in that case, only the integers are sorted 663b2f1ab58SBarry Smith 664b2f1ab58SBarry Smith */ 665b2f1ab58SBarry Smith #undef __FUNCT__ 666b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mergesort" 667b2f1ab58SBarry Smith PetscErrorCode spbas_mergesort(PetscInt nnz, PetscInt *icol, PetscScalar *val) 668b2f1ab58SBarry Smith { 669b2f1ab58SBarry Smith PetscInt istep; // Chunk-sizes of already sorted parts of arrays 670b2f1ab58SBarry Smith PetscInt i, i1, i2; // Loop counters for (partly) sorted arrays 671b2f1ab58SBarry Smith PetscInt istart, i1end, i2end; // start of newly sorted array part, end of both parts 672b2f1ab58SBarry Smith PetscInt *ialloc; // Allocated arrays 673b2f1ab58SBarry Smith PetscScalar *valloc=NULL; 674b2f1ab58SBarry Smith PetscInt *iswap; // auxiliary pointers for swapping 675b2f1ab58SBarry Smith PetscScalar *vswap; 676b2f1ab58SBarry Smith PetscInt *ihlp1; // Pointers to new version of arrays, 677b2f1ab58SBarry Smith PetscScalar *vhlp1=NULL; // (arrays under construction) 678b2f1ab58SBarry Smith PetscInt *ihlp2; // Pointers to previous version of arrays, 679b2f1ab58SBarry Smith PetscScalar *vhlp2=NULL; 680b2f1ab58SBarry Smith PetscErrorCode ierr; 681b2f1ab58SBarry Smith 682b2f1ab58SBarry Smith // Create arrays 683b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&ialloc); CHKERRQ(ierr); 684b2f1ab58SBarry Smith ihlp1 = ialloc; 685b2f1ab58SBarry Smith ihlp2 = icol; 686b2f1ab58SBarry Smith 687b2f1ab58SBarry Smith if (val) 688b2f1ab58SBarry Smith { 689b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscScalar),&valloc); CHKERRQ(ierr); 690b2f1ab58SBarry Smith vhlp1 = valloc; 691b2f1ab58SBarry Smith vhlp2 = val; 692b2f1ab58SBarry Smith } 693b2f1ab58SBarry Smith 694b2f1ab58SBarry Smith 695b2f1ab58SBarry Smith // Sorted array chunks are first 1 long, and increase until they are the 696b2f1ab58SBarry Smith // complete array 697b2f1ab58SBarry Smith for (istep=1; istep<nnz; istep*=2) 698b2f1ab58SBarry Smith { 699b2f1ab58SBarry Smith // 700b2f1ab58SBarry Smith // Combine sorted parts 701b2f1ab58SBarry Smith // istart:istart+istep-1 and istart+istep-1:istart+2*istep-1 702b2f1ab58SBarry Smith // of ihlp2 and vhlp2 703b2f1ab58SBarry Smith // 704b2f1ab58SBarry Smith // into one sorted part 705b2f1ab58SBarry Smith // istart:istart+2*istep-1 706b2f1ab58SBarry Smith // of ihlp1 and vhlp1 707b2f1ab58SBarry Smith // 708b2f1ab58SBarry Smith for (istart=0; istart<nnz; istart+=2*istep) 709b2f1ab58SBarry Smith { 710b2f1ab58SBarry Smith 711b2f1ab58SBarry Smith // Set counters and bound array part endings 712b2f1ab58SBarry Smith i1=istart; i1end = i1+istep; if (i1end>nnz) {i1end=nnz;} 713b2f1ab58SBarry Smith i2=istart+istep; i2end = i2+istep; if (i2end>nnz) {i2end=nnz;} 714b2f1ab58SBarry Smith 715b2f1ab58SBarry Smith // Merge the two array parts 716b2f1ab58SBarry Smith if (val) 717b2f1ab58SBarry Smith { 718b2f1ab58SBarry Smith for (i=istart; i<i2end; i++) 719b2f1ab58SBarry Smith { 720b2f1ab58SBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) 721b2f1ab58SBarry Smith { 722b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 723b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 724b2f1ab58SBarry Smith i1++; 725b2f1ab58SBarry Smith } 726b2f1ab58SBarry Smith else if (i2<i2end ) 727b2f1ab58SBarry Smith { 728b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 729b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i2]; 730b2f1ab58SBarry Smith i2++; 731b2f1ab58SBarry Smith } 732b2f1ab58SBarry Smith else 733b2f1ab58SBarry Smith { 734b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 735b2f1ab58SBarry Smith vhlp1[i] = vhlp2[i1]; 736b2f1ab58SBarry Smith i1++; 737b2f1ab58SBarry Smith } 738b2f1ab58SBarry Smith } 739b2f1ab58SBarry Smith } 740b2f1ab58SBarry Smith else 741b2f1ab58SBarry Smith { 742b2f1ab58SBarry Smith for (i=istart; i<i2end; i++) 743b2f1ab58SBarry Smith { 744b2f1ab58SBarry Smith if (i1<i1end && i2<i2end && ihlp2[i1] < ihlp2[i2]) 745b2f1ab58SBarry Smith { 746b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 747b2f1ab58SBarry Smith i1++; 748b2f1ab58SBarry Smith } 749b2f1ab58SBarry Smith else if (i2<i2end ) 750b2f1ab58SBarry Smith { 751b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i2]; 752b2f1ab58SBarry Smith i2++; 753b2f1ab58SBarry Smith } 754b2f1ab58SBarry Smith else 755b2f1ab58SBarry Smith { 756b2f1ab58SBarry Smith ihlp1[i] = ihlp2[i1]; 757b2f1ab58SBarry Smith i1++; 758b2f1ab58SBarry Smith } 759b2f1ab58SBarry Smith } 760b2f1ab58SBarry Smith } 761b2f1ab58SBarry Smith } 762b2f1ab58SBarry Smith 763b2f1ab58SBarry Smith // Swap the two array sets 764b2f1ab58SBarry Smith iswap = ihlp2; ihlp2 = ihlp1; ihlp1 = iswap; 765b2f1ab58SBarry Smith vswap = vhlp2; vhlp2 = vhlp1; vhlp1 = vswap; 766b2f1ab58SBarry Smith } 767b2f1ab58SBarry Smith 768b2f1ab58SBarry Smith // Copy one more time in case the sorted arrays are the temporary ones 769b2f1ab58SBarry Smith if (ihlp2 != icol) 770b2f1ab58SBarry Smith { 771b2f1ab58SBarry Smith for (i=0; i<nnz; i++) { icol[i] = ihlp2[i]; } 772b2f1ab58SBarry Smith if (val) 773b2f1ab58SBarry Smith { 774b2f1ab58SBarry Smith for (i=0; i<nnz; i++) { val[i] = vhlp2[i]; } 775b2f1ab58SBarry Smith } 776b2f1ab58SBarry Smith } 777b2f1ab58SBarry Smith 778b2f1ab58SBarry Smith // Remove help arrays 779b2f1ab58SBarry Smith ierr = PetscFree(ialloc); CHKERRQ(ierr); 780b2f1ab58SBarry Smith if(val){ierr = PetscFree(valloc); CHKERRQ(ierr);} 781b2f1ab58SBarry Smith PetscFunctionReturn(0); 782b2f1ab58SBarry Smith } 783b2f1ab58SBarry Smith 784b2f1ab58SBarry Smith /* 785b2f1ab58SBarry Smith spbas_apply_reordering_rows: 786b2f1ab58SBarry Smith apply the given reordering to the rows: matrix_A = matrix_A(perm,:); 787b2f1ab58SBarry Smith */ 788b2f1ab58SBarry Smith #undef __FUNCT__ 789b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_rows" 790b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_rows(spbas_matrix *matrix_A, const PetscInt *permutation) 791b2f1ab58SBarry Smith { 792b2f1ab58SBarry Smith PetscInt i,j,ip; 793b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 794b2f1ab58SBarry Smith PetscInt * row_nnz; 795b2f1ab58SBarry Smith PetscInt **icols; 796d36aa34eSBarry Smith PetscTruth do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 797b2f1ab58SBarry Smith PetscScalar **vals=NULL; 798b2f1ab58SBarry Smith PetscErrorCode ierr; 799b2f1ab58SBarry Smith 800b2f1ab58SBarry Smith PetscFunctionBegin; 801b2f1ab58SBarry Smith if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) 802b2f1ab58SBarry Smith { 803b2f1ab58SBarry Smith SETERRQ( PETSC_ERR_SUP_SYS, 804b2f1ab58SBarry Smith "must have diagonal offsets in pattern\n"); 805b2f1ab58SBarry Smith } 806b2f1ab58SBarry Smith 807b2f1ab58SBarry Smith if (do_values) 808b2f1ab58SBarry Smith { 809b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscScalar*)*nrows, &vals); CHKERRQ(ierr); 810b2f1ab58SBarry Smith } 811b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt)*nrows, &row_nnz); CHKERRQ(ierr); 812b2f1ab58SBarry Smith ierr = PetscMalloc( sizeof(PetscInt*)*nrows, &icols); CHKERRQ(ierr); 813b2f1ab58SBarry Smith 814b2f1ab58SBarry Smith for (i=0; i<nrows;i++) 815b2f1ab58SBarry Smith { 816b2f1ab58SBarry Smith ip = permutation[i]; 817b2f1ab58SBarry Smith if (do_values) {vals[i] = matrix_A->values[ip];} 818b2f1ab58SBarry Smith icols[i] = matrix_A->icols[ip]; 819b2f1ab58SBarry Smith row_nnz[i] = matrix_A->row_nnz[ip]; 820b2f1ab58SBarry Smith for (j=0; j<row_nnz[i]; j++) { icols[i][j] += ip-i; } 821b2f1ab58SBarry Smith } 822b2f1ab58SBarry Smith 823b2f1ab58SBarry Smith if (do_values){ ierr = PetscFree(matrix_A->values); CHKERRQ(ierr);} 824b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->icols); CHKERRQ(ierr); 825b2f1ab58SBarry Smith ierr = PetscFree(matrix_A->row_nnz); CHKERRQ(ierr); 826b2f1ab58SBarry Smith 827b2f1ab58SBarry Smith if (do_values) { matrix_A->values = vals; } 828b2f1ab58SBarry Smith matrix_A->icols = icols; 829b2f1ab58SBarry Smith matrix_A->row_nnz = row_nnz; 830b2f1ab58SBarry Smith 831b2f1ab58SBarry Smith PetscFunctionReturn(0); 832b2f1ab58SBarry Smith } 833b2f1ab58SBarry Smith 834b2f1ab58SBarry Smith 835b2f1ab58SBarry Smith /* 836b2f1ab58SBarry Smith spbas_apply_reordering_cols: 837b2f1ab58SBarry Smith apply the given reordering to the columns: matrix_A(:,perm) = matrix_A; 838b2f1ab58SBarry Smith */ 839b2f1ab58SBarry Smith #undef __FUNCT__ 840b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering_cols" 841b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering_cols( spbas_matrix *matrix_A,const PetscInt *permutation) 842b2f1ab58SBarry Smith { 843b2f1ab58SBarry Smith PetscInt i,j; 844b2f1ab58SBarry Smith PetscInt nrows=matrix_A->nrows; 845b2f1ab58SBarry Smith PetscInt row_nnz; 846b2f1ab58SBarry Smith PetscInt *icols; 847d36aa34eSBarry Smith PetscTruth do_values = matrix_A->values ? PETSC_TRUE : PETSC_FALSE; 848b2f1ab58SBarry Smith PetscScalar *vals=NULL; 849b2f1ab58SBarry Smith 850b2f1ab58SBarry Smith PetscErrorCode ierr; 851b2f1ab58SBarry Smith 852b2f1ab58SBarry Smith PetscFunctionBegin; 853b2f1ab58SBarry Smith 854b2f1ab58SBarry Smith if (matrix_A->col_idx_type != SPBAS_DIAGONAL_OFFSETS) 855b2f1ab58SBarry Smith { 856b2f1ab58SBarry Smith SETERRQ( PETSC_ERR_SUP_SYS, 857b2f1ab58SBarry Smith "must have diagonal offsets in pattern\n"); 858b2f1ab58SBarry Smith } 859b2f1ab58SBarry Smith 860b2f1ab58SBarry Smith for (i=0; i<nrows;i++) 861b2f1ab58SBarry Smith { 862b2f1ab58SBarry Smith icols = matrix_A->icols[i]; 863b2f1ab58SBarry Smith row_nnz = matrix_A->row_nnz[i]; 864b2f1ab58SBarry Smith if (do_values) { vals = matrix_A->values[i]; } 865b2f1ab58SBarry Smith 866b2f1ab58SBarry Smith for (j=0; j<row_nnz; j++) 867b2f1ab58SBarry Smith { 868b2f1ab58SBarry Smith icols[j] = permutation[i+icols[j]]-i; 869b2f1ab58SBarry Smith } 870b2f1ab58SBarry Smith ierr = spbas_mergesort(row_nnz, icols, vals); CHKERRQ(ierr); 871b2f1ab58SBarry Smith } 872b2f1ab58SBarry Smith 873b2f1ab58SBarry Smith PetscFunctionReturn(0); 874b2f1ab58SBarry Smith } 875b2f1ab58SBarry Smith 876b2f1ab58SBarry Smith /* 877b2f1ab58SBarry Smith spbas_apply_reordering: 878b2f1ab58SBarry Smith apply the given reordering: matrix_A(perm,perm) = matrix_A; 879b2f1ab58SBarry Smith */ 880b2f1ab58SBarry Smith #undef __FUNCT__ 881b2f1ab58SBarry Smith #define __FUNCT__ "spbas_apply_reordering" 882b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering( spbas_matrix *matrix_A, const PetscInt *permutation, const PetscInt * inv_perm) 883b2f1ab58SBarry Smith { 884b2f1ab58SBarry Smith PetscErrorCode ierr; 885b2f1ab58SBarry Smith PetscFunctionBegin; 886b2f1ab58SBarry Smith ierr = spbas_apply_reordering_rows( matrix_A, inv_perm); CHKERRQ(ierr); 887b2f1ab58SBarry Smith ierr = spbas_apply_reordering_cols( matrix_A, permutation); CHKERRQ(ierr); 888b2f1ab58SBarry Smith PetscFunctionReturn(0); 889b2f1ab58SBarry Smith } 890b2f1ab58SBarry Smith 891b2f1ab58SBarry Smith #undef __FUNCT__ 892b2f1ab58SBarry Smith #define __FUNCT__ "spbas_pattern_only" 893b2f1ab58SBarry Smith PetscErrorCode spbas_pattern_only(PetscInt nrows, PetscInt ncols, PetscInt *ai, PetscInt *aj, spbas_matrix * result) 894b2f1ab58SBarry Smith { 895b2f1ab58SBarry Smith spbas_matrix retval; 896b2f1ab58SBarry Smith PetscInt i, j, i0, r_nnz; 897b2f1ab58SBarry Smith PetscErrorCode ierr; 898b2f1ab58SBarry Smith 899b2f1ab58SBarry Smith PetscFunctionBegin; 900b2f1ab58SBarry Smith 901b2f1ab58SBarry Smith // Copy input values 902b2f1ab58SBarry Smith retval.nrows = nrows; 903b2f1ab58SBarry Smith retval.ncols = ncols; 904b2f1ab58SBarry Smith retval.nnz = ai[nrows]; 905b2f1ab58SBarry Smith 906b2f1ab58SBarry Smith retval.block_data = PETSC_TRUE; 907b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 908b2f1ab58SBarry Smith 909b2f1ab58SBarry Smith // Allocate output matrix 910b2f1ab58SBarry Smith ierr = spbas_allocate_pattern(&retval, PETSC_FALSE); CHKERRQ(ierr); 911b2f1ab58SBarry Smith for (i=0; i<nrows; i++) {retval.row_nnz[i] = ai[i+1]-ai[i];} 912b2f1ab58SBarry Smith ierr = spbas_allocate_data(&retval); CHKERRQ(ierr); 913b2f1ab58SBarry Smith // Copy the structure 914b2f1ab58SBarry Smith for (i = 0; i<retval.nrows; i++) 915b2f1ab58SBarry Smith { 916b2f1ab58SBarry Smith i0 = ai[i]; 917b2f1ab58SBarry Smith r_nnz = ai[i+1]-i0; 918b2f1ab58SBarry Smith 919b2f1ab58SBarry Smith for (j=0; j<r_nnz; j++) 920b2f1ab58SBarry Smith { 921b2f1ab58SBarry Smith retval.icols[i][j] = aj[i0+j]-i; 922b2f1ab58SBarry Smith } 923b2f1ab58SBarry Smith } 924b2f1ab58SBarry Smith 925b2f1ab58SBarry Smith // Set return value 926b2f1ab58SBarry Smith *result = retval; 927b2f1ab58SBarry Smith PetscFunctionReturn(0); 928b2f1ab58SBarry Smith } 929b2f1ab58SBarry Smith 930b2f1ab58SBarry Smith 931b2f1ab58SBarry Smith /* 932b2f1ab58SBarry Smith spbas_mark_row_power: 933b2f1ab58SBarry Smith Mark the columns in row 'row' which are nonzero in 934b2f1ab58SBarry Smith matrix^2log(marker). 935b2f1ab58SBarry Smith */ 936b2f1ab58SBarry Smith #undef __FUNCT__ 937b2f1ab58SBarry Smith #define __FUNCT__ "spbas_mark_row_power" 938b2f1ab58SBarry Smith PetscErrorCode spbas_mark_row_power( 939b2f1ab58SBarry Smith PetscInt *iwork, // marker-vector 940b2f1ab58SBarry Smith PetscInt row, // row for which the columns are marked 941b2f1ab58SBarry Smith spbas_matrix * in_matrix, // matrix for which the power is being 942b2f1ab58SBarry Smith // calculated 943b2f1ab58SBarry Smith PetscInt marker, // marker-value: 2^power 944b2f1ab58SBarry Smith PetscInt minmrk, // lower bound for marked points 945b2f1ab58SBarry Smith PetscInt maxmrk) // upper bound for marked points 946b2f1ab58SBarry Smith { 947b2f1ab58SBarry Smith PetscErrorCode ierr; 948b2f1ab58SBarry Smith PetscInt i,j, nnz; 949b2f1ab58SBarry Smith 950b2f1ab58SBarry Smith PetscFunctionBegin; 951b2f1ab58SBarry Smith nnz = in_matrix->row_nnz[row]; 952b2f1ab58SBarry Smith 953b2f1ab58SBarry Smith // For higher powers, call this function recursively 954b2f1ab58SBarry Smith if (marker>1) 955b2f1ab58SBarry Smith { 956b2f1ab58SBarry Smith for (i=0; i<nnz;i++) 957b2f1ab58SBarry Smith { 958b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 959b2f1ab58SBarry Smith if (minmrk<=j && j<maxmrk && iwork[j] < marker ) 960b2f1ab58SBarry Smith { 961b2f1ab58SBarry Smith ierr = spbas_mark_row_power( iwork, row + in_matrix->icols[row][i], 962b2f1ab58SBarry Smith in_matrix, marker/2,minmrk,maxmrk); 963b2f1ab58SBarry Smith CHKERRQ(ierr); 964b2f1ab58SBarry Smith iwork[j] |= marker; 965b2f1ab58SBarry Smith } 966b2f1ab58SBarry Smith } 967b2f1ab58SBarry Smith } 968b2f1ab58SBarry Smith else 969b2f1ab58SBarry Smith { 970b2f1ab58SBarry Smith // Mark the columns reached 971b2f1ab58SBarry Smith for (i=0; i<nnz;i++) 972b2f1ab58SBarry Smith { 973b2f1ab58SBarry Smith j = row + in_matrix->icols[row][i]; 974b2f1ab58SBarry Smith if (minmrk<=j && j<maxmrk) { iwork[j] |= 1; } 975b2f1ab58SBarry Smith } 976b2f1ab58SBarry Smith } 977b2f1ab58SBarry Smith 978b2f1ab58SBarry Smith PetscFunctionReturn(0); 979b2f1ab58SBarry Smith } 980b2f1ab58SBarry Smith 981b2f1ab58SBarry Smith 982b2f1ab58SBarry Smith /* 983b2f1ab58SBarry Smith spbas_power 984b2f1ab58SBarry Smith Calculate sparseness patterns for incomplete Cholesky decompositions 985b2f1ab58SBarry Smith of a given order: (almost) all nonzeros of the matrix^(order+1) which 986b2f1ab58SBarry Smith are inside the band width are found and stored in the output sparseness 987b2f1ab58SBarry Smith pattern. 988b2f1ab58SBarry Smith */ 989b2f1ab58SBarry Smith #undef __FUNCT__ 990b2f1ab58SBarry Smith #define __FUNCT__ "spbas_power" 991b2f1ab58SBarry Smith PetscErrorCode spbas_power (spbas_matrix in_matrix,PetscInt power, spbas_matrix * result) 992b2f1ab58SBarry Smith { 993b2f1ab58SBarry Smith spbas_matrix retval; 994b2f1ab58SBarry Smith PetscInt nrows = in_matrix.nrows; 995b2f1ab58SBarry Smith PetscInt ncols = in_matrix.ncols; 996b2f1ab58SBarry Smith PetscInt i, j, kend; 997b2f1ab58SBarry Smith PetscInt nnz, inz; 998b2f1ab58SBarry Smith PetscInt *iwork; 999b2f1ab58SBarry Smith PetscInt marker; 1000b2f1ab58SBarry Smith PetscInt maxmrk=0; 1001b2f1ab58SBarry Smith PetscErrorCode ierr; 1002b2f1ab58SBarry Smith 1003b2f1ab58SBarry Smith PetscFunctionBegin; 1004b2f1ab58SBarry Smith if (in_matrix.col_idx_type != SPBAS_DIAGONAL_OFFSETS) 1005b2f1ab58SBarry Smith { 1006b2f1ab58SBarry Smith SETERRQ( PETSC_ERR_SUP_SYS, 1007b2f1ab58SBarry Smith "must have diagonal offsets in pattern\n"); 1008b2f1ab58SBarry Smith } 1009b2f1ab58SBarry Smith 1010b2f1ab58SBarry Smith if (ncols != nrows) 1011b2f1ab58SBarry Smith { 1012b2f1ab58SBarry Smith SETERRQ( PETSC_ERR_ARG_INCOMP, 1013b2f1ab58SBarry Smith "Dimension error\n"); 1014b2f1ab58SBarry Smith } 1015b2f1ab58SBarry Smith 1016b2f1ab58SBarry Smith if (in_matrix.values) 1017b2f1ab58SBarry Smith { 1018b2f1ab58SBarry Smith SETERRQ( PETSC_ERR_ARG_INCOMP, 1019b2f1ab58SBarry Smith "Input array must be sparseness pattern (no values)"); 1020b2f1ab58SBarry Smith } 1021b2f1ab58SBarry Smith 1022b2f1ab58SBarry Smith if (power<=0) 1023b2f1ab58SBarry Smith { 1024b2f1ab58SBarry Smith SETERRQ( PETSC_ERR_SUP_SYS, 1025b2f1ab58SBarry Smith "Power must be 1 or up"); 1026b2f1ab58SBarry Smith } 1027b2f1ab58SBarry Smith 1028b2f1ab58SBarry Smith // Copy input values 1029b2f1ab58SBarry Smith retval.nrows = ncols; 1030b2f1ab58SBarry Smith retval.ncols = nrows; 1031b2f1ab58SBarry Smith retval.nnz = 0; 1032b2f1ab58SBarry Smith retval.col_idx_type = SPBAS_DIAGONAL_OFFSETS; 1033b2f1ab58SBarry Smith retval.block_data = PETSC_FALSE; 1034b2f1ab58SBarry Smith 1035b2f1ab58SBarry Smith // Allocate sparseness pattern 1036d36aa34eSBarry Smith ierr = spbas_allocate_pattern(&retval, in_matrix.values ? PETSC_TRUE : PETSC_FALSE); 1037b2f1ab58SBarry Smith CHKERRQ(ierr); 1038b2f1ab58SBarry Smith 1039b2f1ab58SBarry Smith // Allocate marker array 1040b2f1ab58SBarry Smith ierr = PetscMalloc(nrows * sizeof(PetscInt), &iwork); CHKERRQ(ierr); 1041b2f1ab58SBarry Smith 1042b2f1ab58SBarry Smith // Erase the pattern for this row 1043b2f1ab58SBarry Smith memset( (void *) iwork, 0, retval.nrows*sizeof(PetscInt)); 1044b2f1ab58SBarry Smith 1045b2f1ab58SBarry Smith // Calculate marker values 1046b2f1ab58SBarry Smith marker = 1; for (i=1; i<power; i++) {marker*=2;} 1047b2f1ab58SBarry Smith 1048b2f1ab58SBarry Smith for (i=0; i<nrows; i++) 1049b2f1ab58SBarry Smith { 1050b2f1ab58SBarry Smith // Calculate the pattern for each row 1051b2f1ab58SBarry Smith 1052b2f1ab58SBarry Smith nnz = in_matrix.row_nnz[i]; 1053b2f1ab58SBarry Smith kend = i+in_matrix.icols[i][nnz-1]; 1054b2f1ab58SBarry Smith if (maxmrk<=kend) {maxmrk=kend+1;} 1055b2f1ab58SBarry Smith ierr = spbas_mark_row_power( iwork, i, &in_matrix, marker, 1056b2f1ab58SBarry Smith i, maxmrk); CHKERRQ(ierr); 1057b2f1ab58SBarry Smith 1058b2f1ab58SBarry Smith // Count the columns 1059b2f1ab58SBarry Smith nnz = 0; 1060b2f1ab58SBarry Smith for (j=i; j<maxmrk; j++) {nnz+= (iwork[j]!=0);} 1061b2f1ab58SBarry Smith 1062b2f1ab58SBarry Smith // Allocate the column indices 1063b2f1ab58SBarry Smith retval.row_nnz[i] = nnz; 1064b2f1ab58SBarry Smith ierr = PetscMalloc(nnz*sizeof(PetscInt),&retval.icols[i]); CHKERRQ(ierr); 1065b2f1ab58SBarry Smith 1066b2f1ab58SBarry Smith // Administrate the column indices 1067b2f1ab58SBarry Smith inz = 0; 1068b2f1ab58SBarry Smith for (j=i; j<maxmrk; j++) 1069b2f1ab58SBarry Smith { 1070b2f1ab58SBarry Smith if (iwork[j]) 1071b2f1ab58SBarry Smith { 1072b2f1ab58SBarry Smith retval.icols[i][inz] = j-i; 1073b2f1ab58SBarry Smith inz++; 1074b2f1ab58SBarry Smith iwork[j]=0; 1075b2f1ab58SBarry Smith } 1076b2f1ab58SBarry Smith } 1077b2f1ab58SBarry Smith retval.nnz += nnz; 1078b2f1ab58SBarry Smith }; 1079b2f1ab58SBarry Smith 1080b2f1ab58SBarry Smith ierr = PetscFree(iwork); CHKERRQ(ierr); 1081b2f1ab58SBarry Smith 1082b2f1ab58SBarry Smith *result = retval; 1083b2f1ab58SBarry Smith PetscFunctionReturn(0); 1084b2f1ab58SBarry Smith } 1085b2f1ab58SBarry Smith 1086b2f1ab58SBarry Smith 1087b2f1ab58SBarry Smith 1088b2f1ab58SBarry Smith /* 1089b2f1ab58SBarry Smith spbas_keep_upper: 1090b2f1ab58SBarry Smith remove the lower part of the matrix: keep the upper part 1091b2f1ab58SBarry Smith */ 1092b2f1ab58SBarry Smith #undef __FUNCT__ 1093b2f1ab58SBarry Smith #define __FUNCT__ "spbas_keep_upper" 1094b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper( spbas_matrix * inout_matrix) 1095b2f1ab58SBarry Smith { 1096b2f1ab58SBarry Smith PetscInt i, j; 1097b2f1ab58SBarry Smith PetscInt jstart; 1098b2f1ab58SBarry Smith 1099b2f1ab58SBarry Smith PetscFunctionBegin; 1100b2f1ab58SBarry Smith if (inout_matrix->block_data) 1101b2f1ab58SBarry Smith { 1102b2f1ab58SBarry Smith SETERRQ( PETSC_ERR_SUP_SYS, 1103b2f1ab58SBarry Smith "Not yet for block data matrices\n"); 1104b2f1ab58SBarry Smith } 1105b2f1ab58SBarry Smith 1106b2f1ab58SBarry Smith for (i=0; i<inout_matrix->nrows; i++) 1107b2f1ab58SBarry Smith { 1108b2f1ab58SBarry Smith for (jstart=0; 1109b2f1ab58SBarry Smith (jstart<inout_matrix->row_nnz[i]) && 1110b2f1ab58SBarry Smith (inout_matrix->icols[i][jstart]<0); jstart++){} 1111b2f1ab58SBarry Smith 1112b2f1ab58SBarry Smith if (jstart>0) 1113b2f1ab58SBarry Smith { 1114b2f1ab58SBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) 1115b2f1ab58SBarry Smith { 1116b2f1ab58SBarry Smith inout_matrix->icols[i][j] = 1117b2f1ab58SBarry Smith inout_matrix->icols[i][j+jstart]; 1118b2f1ab58SBarry Smith } 1119b2f1ab58SBarry Smith 1120b2f1ab58SBarry Smith if (inout_matrix->values) 1121b2f1ab58SBarry Smith { 1122b2f1ab58SBarry Smith for (j=0; j<inout_matrix->row_nnz[i]-jstart; j++) 1123b2f1ab58SBarry Smith { 1124b2f1ab58SBarry Smith inout_matrix->values[i][j] = 1125b2f1ab58SBarry Smith inout_matrix->values[i][j+jstart]; 1126b2f1ab58SBarry Smith } 1127b2f1ab58SBarry Smith } 1128b2f1ab58SBarry Smith 1129b2f1ab58SBarry Smith inout_matrix->row_nnz[i] -= jstart; 1130b2f1ab58SBarry Smith 1131b2f1ab58SBarry Smith inout_matrix->icols[i] = (PetscInt*) realloc( 1132b2f1ab58SBarry Smith (void*) inout_matrix->icols[i], 1133b2f1ab58SBarry Smith inout_matrix->row_nnz[i]*sizeof(PetscInt)); 1134b2f1ab58SBarry Smith 1135b2f1ab58SBarry Smith if (inout_matrix->values) 1136b2f1ab58SBarry Smith { 1137b2f1ab58SBarry Smith inout_matrix->values[i] = (PetscScalar*) realloc( 1138b2f1ab58SBarry Smith (void*) inout_matrix->values[i], 1139b2f1ab58SBarry Smith inout_matrix->row_nnz[i]*sizeof(PetscScalar)); 1140b2f1ab58SBarry Smith } 1141b2f1ab58SBarry Smith inout_matrix->nnz -= jstart; 1142b2f1ab58SBarry Smith } 1143b2f1ab58SBarry Smith } 1144b2f1ab58SBarry Smith PetscFunctionReturn(0); 1145b2f1ab58SBarry Smith } 1146b2f1ab58SBarry Smith 1147