1*49b5e25fSSatish Balay /*$Id: sbaij.c,v 1.1 2000/01/18 16:45:51 hzhang Exp hzhang $*/ 2*49b5e25fSSatish Balay 3*49b5e25fSSatish Balay /* 4*49b5e25fSSatish Balay Defines the basic matrix operations for the BAIJ (compressed row) 5*49b5e25fSSatish Balay matrix storage format. 6*49b5e25fSSatish Balay */ 7*49b5e25fSSatish Balay #include "petscsys.h" 8*49b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h" 9*49b5e25fSSatish Balay #include "src/vec/vecimpl.h" 10*49b5e25fSSatish Balay #include "src/inline/spops.h" 11*49b5e25fSSatish Balay #include "sbaij.h" 12*49b5e25fSSatish Balay 13*49b5e25fSSatish Balay #define CHUNKSIZE 10 14*49b5e25fSSatish Balay 15*49b5e25fSSatish Balay /* 16*49b5e25fSSatish Balay Checks for missing diagonals 17*49b5e25fSSatish Balay */ 18*49b5e25fSSatish Balay #undef __FUNC__ 19*49b5e25fSSatish Balay #define __FUNC__ "MatMissingDiagonal_SeqSBAIJ" 20*49b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A) 21*49b5e25fSSatish Balay { 22*49b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 23*49b5e25fSSatish Balay int *diag,*jj = a->j,i,ierr; 24*49b5e25fSSatish Balay 25*49b5e25fSSatish Balay PetscFunctionBegin; 26*49b5e25fSSatish Balay ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 27*49b5e25fSSatish Balay diag = a->diag; 28*49b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 29*49b5e25fSSatish Balay if (jj[diag[i]] != i) { 30*49b5e25fSSatish Balay SETERRQ1(1,1,"Matrix is missing diagonal number %d",i); 31*49b5e25fSSatish Balay } 32*49b5e25fSSatish Balay } 33*49b5e25fSSatish Balay PetscFunctionReturn(0); 34*49b5e25fSSatish Balay } 35*49b5e25fSSatish Balay 36*49b5e25fSSatish Balay #undef __FUNC__ 37*49b5e25fSSatish Balay #define __FUNC__ "MatMarkDiagonal_SeqSBAIJ" 38*49b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A) 39*49b5e25fSSatish Balay { 40*49b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 41*49b5e25fSSatish Balay int i,j,*diag,m = a->mbs; 42*49b5e25fSSatish Balay 43*49b5e25fSSatish Balay PetscFunctionBegin; 44*49b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 45*49b5e25fSSatish Balay 46*49b5e25fSSatish Balay diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag); 47*49b5e25fSSatish Balay PLogObjectMemory(A,(m+1)*sizeof(int)); 48*49b5e25fSSatish Balay for (i=0; i<m; i++) { 49*49b5e25fSSatish Balay diag[i] = a->i[i+1]; 50*49b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 51*49b5e25fSSatish Balay if (a->j[j] == i) { 52*49b5e25fSSatish Balay diag[i] = j; 53*49b5e25fSSatish Balay break; 54*49b5e25fSSatish Balay } 55*49b5e25fSSatish Balay } 56*49b5e25fSSatish Balay } 57*49b5e25fSSatish Balay a->diag = diag; 58*49b5e25fSSatish Balay PetscFunctionReturn(0); 59*49b5e25fSSatish Balay } 60*49b5e25fSSatish Balay 61*49b5e25fSSatish Balay 62*49b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 63*49b5e25fSSatish Balay 64*49b5e25fSSatish Balay #undef __FUNC__ 65*49b5e25fSSatish Balay #define __FUNC__ "MatGetRowIJ_SeqSBAIJ" 66*49b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 67*49b5e25fSSatish Balay { 68*49b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 69*49b5e25fSSatish Balay int ierr,n = a->mbs,i; 70*49b5e25fSSatish Balay 71*49b5e25fSSatish Balay PetscFunctionBegin; 72*49b5e25fSSatish Balay *nn = n; 73*49b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 74*49b5e25fSSatish Balay if (symmetric) { 75*49b5e25fSSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 76*49b5e25fSSatish Balay } else if (oshift == 1) { 77*49b5e25fSSatish Balay /* temporarily add 1 to i and j indices */ 78*49b5e25fSSatish Balay int nz = a->i[n] + 1; 79*49b5e25fSSatish Balay for (i=0; i<nz; i++) a->j[i]++; 80*49b5e25fSSatish Balay for (i=0; i<n+1; i++) a->i[i]++; 81*49b5e25fSSatish Balay *ia = a->i; *ja = a->j; 82*49b5e25fSSatish Balay } else { 83*49b5e25fSSatish Balay *ia = a->i; *ja = a->j; 84*49b5e25fSSatish Balay } 85*49b5e25fSSatish Balay 86*49b5e25fSSatish Balay PetscFunctionReturn(0); 87*49b5e25fSSatish Balay } 88*49b5e25fSSatish Balay 89*49b5e25fSSatish Balay #undef __FUNC__ 90*49b5e25fSSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqSBAIJ" 91*49b5e25fSSatish Balay static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 92*49b5e25fSSatish Balay PetscTruth *done) 93*49b5e25fSSatish Balay { 94*49b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 95*49b5e25fSSatish Balay int i,n = a->mbs,ierr; 96*49b5e25fSSatish Balay 97*49b5e25fSSatish Balay PetscFunctionBegin; 98*49b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 99*49b5e25fSSatish Balay if (symmetric) { 100*49b5e25fSSatish Balay ierr = PetscFree(*ia);CHKERRQ(ierr); 101*49b5e25fSSatish Balay ierr = PetscFree(*ja);CHKERRQ(ierr); 102*49b5e25fSSatish Balay } else if (oshift == 1) { 103*49b5e25fSSatish Balay int nz = a->i[n]; 104*49b5e25fSSatish Balay for (i=0; i<nz; i++) a->j[i]--; 105*49b5e25fSSatish Balay for (i=0; i<n+1; i++) a->i[i]--; 106*49b5e25fSSatish Balay } 107*49b5e25fSSatish Balay PetscFunctionReturn(0); 108*49b5e25fSSatish Balay } 109*49b5e25fSSatish Balay 110*49b5e25fSSatish Balay #undef __FUNC__ 111*49b5e25fSSatish Balay #define __FUNC__ "MatGetBlockSize_SeqSBAIJ" 112*49b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 113*49b5e25fSSatish Balay { 114*49b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 115*49b5e25fSSatish Balay 116*49b5e25fSSatish Balay PetscFunctionBegin; 117*49b5e25fSSatish Balay *bs = sbaij->bs; 118*49b5e25fSSatish Balay PetscFunctionReturn(0); 119*49b5e25fSSatish Balay } 120*49b5e25fSSatish Balay 121*49b5e25fSSatish Balay #undef __FUNC__ 122*49b5e25fSSatish Balay #define __FUNC__ "MatDestroy_SeqSBAIJ" 123*49b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A) 124*49b5e25fSSatish Balay { 125*49b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 126*49b5e25fSSatish Balay int ierr; 127*49b5e25fSSatish Balay 128*49b5e25fSSatish Balay PetscFunctionBegin; 129*49b5e25fSSatish Balay if (A->mapping) { 130*49b5e25fSSatish Balay ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); 131*49b5e25fSSatish Balay } 132*49b5e25fSSatish Balay if (A->bmapping) { 133*49b5e25fSSatish Balay ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); 134*49b5e25fSSatish Balay } 135*49b5e25fSSatish Balay if (A->rmap) { 136*49b5e25fSSatish Balay ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 137*49b5e25fSSatish Balay } 138*49b5e25fSSatish Balay if (A->cmap) { 139*49b5e25fSSatish Balay ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 140*49b5e25fSSatish Balay } 141*49b5e25fSSatish Balay #if defined(PETSC_USE_LOG) 142*49b5e25fSSatish Balay PLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",a->m,a->s_nz); 143*49b5e25fSSatish Balay #endif 144*49b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 145*49b5e25fSSatish Balay if (!a->singlemalloc) { 146*49b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 147*49b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 148*49b5e25fSSatish Balay } 149*49b5e25fSSatish Balay if (a->row) { 150*49b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 151*49b5e25fSSatish Balay } 152*49b5e25fSSatish Balay if (a->col) { 153*49b5e25fSSatish Balay ierr = ISDestroy(a->col);CHKERRQ(ierr); 154*49b5e25fSSatish Balay } 155*49b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 156*49b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 157*49b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 158*49b5e25fSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 159*49b5e25fSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 160*49b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 161*49b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 162*49b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 163*49b5e25fSSatish Balay PLogObjectDestroy(A); 164*49b5e25fSSatish Balay PetscHeaderDestroy(A); 165*49b5e25fSSatish Balay PetscFunctionReturn(0); 166*49b5e25fSSatish Balay } 167*49b5e25fSSatish Balay 168*49b5e25fSSatish Balay #undef __FUNC__ 169*49b5e25fSSatish Balay #define __FUNC__ "MatSetOption_SeqSBAIJ" 170*49b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op) 171*49b5e25fSSatish Balay { 172*49b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 173*49b5e25fSSatish Balay 174*49b5e25fSSatish Balay PetscFunctionBegin; 175*49b5e25fSSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = PETSC_TRUE; 176*49b5e25fSSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = PETSC_FALSE; 177*49b5e25fSSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = PETSC_TRUE; 178*49b5e25fSSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = PETSC_FALSE; 179*49b5e25fSSatish Balay else if (op == MAT_KEEP_ZEROED_ROWS) a->keepzeroedrows = PETSC_TRUE; 180*49b5e25fSSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 181*49b5e25fSSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 182*49b5e25fSSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 183*49b5e25fSSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 184*49b5e25fSSatish Balay else if (op == MAT_ROWS_SORTED || 185*49b5e25fSSatish Balay op == MAT_ROWS_UNSORTED || 186*49b5e25fSSatish Balay op == MAT_SYMMETRIC || 187*49b5e25fSSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 188*49b5e25fSSatish Balay op == MAT_YES_NEW_DIAGONALS || 189*49b5e25fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 190*49b5e25fSSatish Balay op == MAT_USE_HASH_TABLE) { 191*49b5e25fSSatish Balay PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 192*49b5e25fSSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 193*49b5e25fSSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 194*49b5e25fSSatish Balay } else { 195*49b5e25fSSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 196*49b5e25fSSatish Balay } 197*49b5e25fSSatish Balay PetscFunctionReturn(0); 198*49b5e25fSSatish Balay } 199*49b5e25fSSatish Balay 200*49b5e25fSSatish Balay 201*49b5e25fSSatish Balay #undef __FUNC__ 202*49b5e25fSSatish Balay #define __FUNC__ "MatGetSize_SeqSBAIJ" 203*49b5e25fSSatish Balay int MatGetSize_SeqSBAIJ(Mat A,int *m,int *n) 204*49b5e25fSSatish Balay { 205*49b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 206*49b5e25fSSatish Balay 207*49b5e25fSSatish Balay PetscFunctionBegin; 208*49b5e25fSSatish Balay if (m) *m = a->m; 209*49b5e25fSSatish Balay if (n) *n = a->m; 210*49b5e25fSSatish Balay PetscFunctionReturn(0); 211*49b5e25fSSatish Balay } 212*49b5e25fSSatish Balay 213*49b5e25fSSatish Balay #undef __FUNC__ 214*49b5e25fSSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqSBAIJ" 215*49b5e25fSSatish Balay int MatGetOwnershipRange_SeqSBAIJ(Mat A,int *m,int *n) 216*49b5e25fSSatish Balay { 217*49b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 218*49b5e25fSSatish Balay 219*49b5e25fSSatish Balay PetscFunctionBegin; 220*49b5e25fSSatish Balay *m = 0; *n = a->m; 221*49b5e25fSSatish Balay PetscFunctionReturn(0); 222*49b5e25fSSatish Balay } 223*49b5e25fSSatish Balay 224*49b5e25fSSatish Balay #undef __FUNC__ 225*49b5e25fSSatish Balay #define __FUNC__ "MatGetRow_SeqSBAIJ" 226*49b5e25fSSatish Balay int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,Scalar **v) 227*49b5e25fSSatish Balay { 228*49b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 229*49b5e25fSSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 230*49b5e25fSSatish Balay MatScalar *aa,*aa_i; 231*49b5e25fSSatish Balay Scalar *v_i; 232*49b5e25fSSatish Balay 233*49b5e25fSSatish Balay PetscFunctionBegin; 234*49b5e25fSSatish Balay bs = a->bs; 235*49b5e25fSSatish Balay ai = a->i; 236*49b5e25fSSatish Balay aj = a->j; 237*49b5e25fSSatish Balay aa = a->a; 238*49b5e25fSSatish Balay bs2 = a->bs2; 239*49b5e25fSSatish Balay 240*49b5e25fSSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 241*49b5e25fSSatish Balay 242*49b5e25fSSatish Balay bn = row/bs; /* Block number */ 243*49b5e25fSSatish Balay bp = row % bs; /* Block position */ 244*49b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 245*49b5e25fSSatish Balay *ncols = bs*M; 246*49b5e25fSSatish Balay 247*49b5e25fSSatish Balay if (v) { 248*49b5e25fSSatish Balay *v = 0; 249*49b5e25fSSatish Balay if (*ncols) { 250*49b5e25fSSatish Balay *v = (Scalar*)PetscMalloc((*ncols+row)*sizeof(Scalar));CHKPTRQ(*v); 251*49b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 252*49b5e25fSSatish Balay v_i = *v + i*bs; 253*49b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 254*49b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 255*49b5e25fSSatish Balay } 256*49b5e25fSSatish Balay } 257*49b5e25fSSatish Balay } 258*49b5e25fSSatish Balay 259*49b5e25fSSatish Balay if (cols) { 260*49b5e25fSSatish Balay *cols = 0; 261*49b5e25fSSatish Balay if (*ncols) { 262*49b5e25fSSatish Balay *cols = (int*)PetscMalloc((*ncols+row)*sizeof(int));CHKPTRQ(*cols); 263*49b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 264*49b5e25fSSatish Balay cols_i = *cols + i*bs; 265*49b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 266*49b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 267*49b5e25fSSatish Balay } 268*49b5e25fSSatish Balay } 269*49b5e25fSSatish Balay } 270*49b5e25fSSatish Balay 271*49b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 272*49b5e25fSSatish Balay v_i = *v + M*bs; 273*49b5e25fSSatish Balay cols_i = *cols + M*bs; 274*49b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 275*49b5e25fSSatish Balay M = ai[i+1] - ai[i]; 276*49b5e25fSSatish Balay for (j=0; j<M; j++){ 277*49b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 278*49b5e25fSSatish Balay if (itmp == bn){ 279*49b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 280*49b5e25fSSatish Balay for (k=0; k<bs; k++) { 281*49b5e25fSSatish Balay *cols_i++ = i*bs+k; 282*49b5e25fSSatish Balay *v_i++ = aa_i[k]; 283*49b5e25fSSatish Balay } 284*49b5e25fSSatish Balay *ncols += bs; 285*49b5e25fSSatish Balay break; 286*49b5e25fSSatish Balay } 287*49b5e25fSSatish Balay } 288*49b5e25fSSatish Balay } 289*49b5e25fSSatish Balay 290*49b5e25fSSatish Balay PetscFunctionReturn(0); 291*49b5e25fSSatish Balay } 292*49b5e25fSSatish Balay 293*49b5e25fSSatish Balay #undef __FUNC__ 294*49b5e25fSSatish Balay #define __FUNC__ "MatRestoreRow_SeqSBAIJ" 295*49b5e25fSSatish Balay int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 296*49b5e25fSSatish Balay { 297*49b5e25fSSatish Balay int ierr; 298*49b5e25fSSatish Balay 299*49b5e25fSSatish Balay PetscFunctionBegin; 300*49b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 301*49b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 302*49b5e25fSSatish Balay PetscFunctionReturn(0); 303*49b5e25fSSatish Balay } 304*49b5e25fSSatish Balay 305*49b5e25fSSatish Balay #undef __FUNC__ 306*49b5e25fSSatish Balay #define __FUNC__ "MatTranspose_SeqSBAIJ" 307*49b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B) 308*49b5e25fSSatish Balay { 309*49b5e25fSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 310*49b5e25fSSatish Balay Mat C; 311*49b5e25fSSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 312*49b5e25fSSatish Balay int *rows,*cols,bs2=a->bs2,refct; 313*49b5e25fSSatish Balay MatScalar *array = a->a; 314*49b5e25fSSatish Balay 315*49b5e25fSSatish Balay PetscFunctionBegin; 316*49b5e25fSSatish Balay if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 317*49b5e25fSSatish Balay col = (int*)PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col); 318*49b5e25fSSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 319*49b5e25fSSatish Balay 320*49b5e25fSSatish Balay for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1; 321*49b5e25fSSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 322*49b5e25fSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 323*49b5e25fSSatish Balay rows = (int*)PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows); 324*49b5e25fSSatish Balay cols = rows + bs; 325*49b5e25fSSatish Balay for (i=0; i<mbs; i++) { 326*49b5e25fSSatish Balay cols[0] = i*bs; 327*49b5e25fSSatish Balay for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1; 328*49b5e25fSSatish Balay len = ai[i+1] - ai[i]; 329*49b5e25fSSatish Balay for (j=0; j<len; j++) { 330*49b5e25fSSatish Balay rows[0] = (*aj++)*bs; 331*49b5e25fSSatish Balay for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1; 332*49b5e25fSSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 333*49b5e25fSSatish Balay array += bs2; 334*49b5e25fSSatish Balay } 335*49b5e25fSSatish Balay } 336*49b5e25fSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 337*49b5e25fSSatish Balay 338*49b5e25fSSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 339*49b5e25fSSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 340*49b5e25fSSatish Balay 341*49b5e25fSSatish Balay if (B) { 342*49b5e25fSSatish Balay *B = C; 343*49b5e25fSSatish Balay } else { 344*49b5e25fSSatish Balay PetscOps *Abops; 345*49b5e25fSSatish Balay MatOps Aops; 346*49b5e25fSSatish Balay 347*49b5e25fSSatish Balay /* This isn't really an in-place transpose */ 348*49b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 349*49b5e25fSSatish Balay if (!a->singlemalloc) { 350*49b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 351*49b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 352*49b5e25fSSatish Balay } 353*49b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 354*49b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 355*49b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 356*49b5e25fSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 357*49b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 358*49b5e25fSSatish Balay 359*49b5e25fSSatish Balay 360*49b5e25fSSatish Balay ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 361*49b5e25fSSatish Balay ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 362*49b5e25fSSatish Balay 363*49b5e25fSSatish Balay /* 364*49b5e25fSSatish Balay This is horrible, horrible code. We need to keep the 365*49b5e25fSSatish Balay A pointers for the bops and ops but copy everything 366*49b5e25fSSatish Balay else from C. 367*49b5e25fSSatish Balay */ 368*49b5e25fSSatish Balay Abops = A->bops; 369*49b5e25fSSatish Balay Aops = A->ops; 370*49b5e25fSSatish Balay refct = A->refct; 371*49b5e25fSSatish Balay ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 372*49b5e25fSSatish Balay A->bops = Abops; 373*49b5e25fSSatish Balay A->ops = Aops; 374*49b5e25fSSatish Balay A->qlist = 0; 375*49b5e25fSSatish Balay A->refct = refct; 376*49b5e25fSSatish Balay 377*49b5e25fSSatish Balay PetscHeaderDestroy(C); 378*49b5e25fSSatish Balay } 379*49b5e25fSSatish Balay PetscFunctionReturn(0); 380*49b5e25fSSatish Balay } 381*49b5e25fSSatish Balay 382*49b5e25fSSatish Balay #undef __FUNC__ 383*49b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Binary" 384*49b5e25fSSatish Balay static int MatView_SeqSBAIJ_Binary(Mat A,Viewer viewer) 385*49b5e25fSSatish Balay { 386*49b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 387*49b5e25fSSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 388*49b5e25fSSatish Balay Scalar *aa; 389*49b5e25fSSatish Balay FILE *file; 390*49b5e25fSSatish Balay 391*49b5e25fSSatish Balay PetscFunctionBegin; 392*49b5e25fSSatish Balay ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 393*49b5e25fSSatish Balay col_lens = (int*)PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 394*49b5e25fSSatish Balay col_lens[0] = MAT_COOKIE; 395*49b5e25fSSatish Balay 396*49b5e25fSSatish Balay col_lens[1] = a->m; 397*49b5e25fSSatish Balay col_lens[2] = a->m; 398*49b5e25fSSatish Balay col_lens[3] = a->s_nz*bs2; 399*49b5e25fSSatish Balay 400*49b5e25fSSatish Balay /* store lengths of each row and write (including header) to file */ 401*49b5e25fSSatish Balay count = 0; 402*49b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 403*49b5e25fSSatish Balay for (j=0; j<bs; j++) { 404*49b5e25fSSatish Balay col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 405*49b5e25fSSatish Balay } 406*49b5e25fSSatish Balay } 407*49b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr); 408*49b5e25fSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 409*49b5e25fSSatish Balay 410*49b5e25fSSatish Balay /* store column indices (zero start index) */ 411*49b5e25fSSatish Balay jj = (int*)PetscMalloc((a->s_nz+1)*bs2*sizeof(int));CHKPTRQ(jj); 412*49b5e25fSSatish Balay count = 0; 413*49b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 414*49b5e25fSSatish Balay for (j=0; j<bs; j++) { 415*49b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 416*49b5e25fSSatish Balay for (l=0; l<bs; l++) { 417*49b5e25fSSatish Balay jj[count++] = bs*a->j[k] + l; 418*49b5e25fSSatish Balay } 419*49b5e25fSSatish Balay } 420*49b5e25fSSatish Balay } 421*49b5e25fSSatish Balay } 422*49b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr); 423*49b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 424*49b5e25fSSatish Balay 425*49b5e25fSSatish Balay /* store nonzero values */ 426*49b5e25fSSatish Balay aa = (Scalar*)PetscMalloc((a->s_nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 427*49b5e25fSSatish Balay count = 0; 428*49b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 429*49b5e25fSSatish Balay for (j=0; j<bs; j++) { 430*49b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 431*49b5e25fSSatish Balay for (l=0; l<bs; l++) { 432*49b5e25fSSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 433*49b5e25fSSatish Balay } 434*49b5e25fSSatish Balay } 435*49b5e25fSSatish Balay } 436*49b5e25fSSatish Balay } 437*49b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr); 438*49b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 439*49b5e25fSSatish Balay 440*49b5e25fSSatish Balay ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 441*49b5e25fSSatish Balay if (file) { 442*49b5e25fSSatish Balay fprintf(file,"-matload_block_size %d\n",a->bs); 443*49b5e25fSSatish Balay } 444*49b5e25fSSatish Balay PetscFunctionReturn(0); 445*49b5e25fSSatish Balay } 446*49b5e25fSSatish Balay 447*49b5e25fSSatish Balay #undef __FUNC__ 448*49b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_ASCII" 449*49b5e25fSSatish Balay static int MatView_SeqSBAIJ_ASCII(Mat A,Viewer viewer) 450*49b5e25fSSatish Balay { 451*49b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 452*49b5e25fSSatish Balay int ierr,i,j,format,bs = a->bs,k,l,bs2=a->bs2; 453*49b5e25fSSatish Balay char *outputname; 454*49b5e25fSSatish Balay 455*49b5e25fSSatish Balay PetscFunctionBegin; 456*49b5e25fSSatish Balay ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 457*49b5e25fSSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 458*49b5e25fSSatish Balay if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 459*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 460*49b5e25fSSatish Balay } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 461*49b5e25fSSatish Balay SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 462*49b5e25fSSatish Balay } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 463*49b5e25fSSatish Balay ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 464*49b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 465*49b5e25fSSatish Balay for (j=0; j<bs; j++) { 466*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 467*49b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 468*49b5e25fSSatish Balay for (l=0; l<bs; l++) { 469*49b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 470*49b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 471*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 472*49b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 473*49b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 474*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 475*49b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 476*49b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 477*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 478*49b5e25fSSatish Balay } 479*49b5e25fSSatish Balay #else 480*49b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 481*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 482*49b5e25fSSatish Balay } 483*49b5e25fSSatish Balay #endif 484*49b5e25fSSatish Balay } 485*49b5e25fSSatish Balay } 486*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 487*49b5e25fSSatish Balay } 488*49b5e25fSSatish Balay } 489*49b5e25fSSatish Balay ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 490*49b5e25fSSatish Balay } else { 491*49b5e25fSSatish Balay ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 492*49b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 493*49b5e25fSSatish Balay for (j=0; j<bs; j++) { 494*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 495*49b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 496*49b5e25fSSatish Balay for (l=0; l<bs; l++) { 497*49b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 498*49b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 499*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 500*49b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 501*49b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 502*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 503*49b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 504*49b5e25fSSatish Balay } else { 505*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 506*49b5e25fSSatish Balay } 507*49b5e25fSSatish Balay #else 508*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 509*49b5e25fSSatish Balay #endif 510*49b5e25fSSatish Balay } 511*49b5e25fSSatish Balay } 512*49b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 513*49b5e25fSSatish Balay } 514*49b5e25fSSatish Balay } 515*49b5e25fSSatish Balay ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 516*49b5e25fSSatish Balay } 517*49b5e25fSSatish Balay ierr = ViewerFlush(viewer);CHKERRQ(ierr); 518*49b5e25fSSatish Balay PetscFunctionReturn(0); 519*49b5e25fSSatish Balay } 520*49b5e25fSSatish Balay 521*49b5e25fSSatish Balay #undef __FUNC__ 522*49b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw_Zoom" 523*49b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw_Zoom(Draw draw,void *Aa) 524*49b5e25fSSatish Balay { 525*49b5e25fSSatish Balay Mat A = (Mat) Aa; 526*49b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 527*49b5e25fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 528*49b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 529*49b5e25fSSatish Balay MatScalar *aa; 530*49b5e25fSSatish Balay MPI_Comm comm; 531*49b5e25fSSatish Balay Viewer viewer; 532*49b5e25fSSatish Balay 533*49b5e25fSSatish Balay PetscFunctionBegin; 534*49b5e25fSSatish Balay /* 535*49b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 536*49b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 537*49b5e25fSSatish Balay rest should return immediately. 538*49b5e25fSSatish Balay */ 539*49b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 540*49b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 541*49b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 542*49b5e25fSSatish Balay 543*49b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 544*49b5e25fSSatish Balay 545*49b5e25fSSatish Balay ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 546*49b5e25fSSatish Balay DrawString(draw, .3*(xl+xr), .3*(yl+yr), DRAW_BLACK, "symmetric"); 547*49b5e25fSSatish Balay 548*49b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 549*49b5e25fSSatish Balay color = DRAW_BLUE; 550*49b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 551*49b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 552*49b5e25fSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 553*49b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 554*49b5e25fSSatish Balay aa = a->a + j*bs2; 555*49b5e25fSSatish Balay for (k=0; k<bs; k++) { 556*49b5e25fSSatish Balay for (l=0; l<bs; l++) { 557*49b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 558*49b5e25fSSatish Balay ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 559*49b5e25fSSatish Balay } 560*49b5e25fSSatish Balay } 561*49b5e25fSSatish Balay } 562*49b5e25fSSatish Balay } 563*49b5e25fSSatish Balay color = DRAW_CYAN; 564*49b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 565*49b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 566*49b5e25fSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 567*49b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 568*49b5e25fSSatish Balay aa = a->a + j*bs2; 569*49b5e25fSSatish Balay for (k=0; k<bs; k++) { 570*49b5e25fSSatish Balay for (l=0; l<bs; l++) { 571*49b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 572*49b5e25fSSatish Balay ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 573*49b5e25fSSatish Balay } 574*49b5e25fSSatish Balay } 575*49b5e25fSSatish Balay } 576*49b5e25fSSatish Balay } 577*49b5e25fSSatish Balay 578*49b5e25fSSatish Balay color = DRAW_RED; 579*49b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 580*49b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 581*49b5e25fSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 582*49b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 583*49b5e25fSSatish Balay aa = a->a + j*bs2; 584*49b5e25fSSatish Balay for (k=0; k<bs; k++) { 585*49b5e25fSSatish Balay for (l=0; l<bs; l++) { 586*49b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 587*49b5e25fSSatish Balay ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 588*49b5e25fSSatish Balay } 589*49b5e25fSSatish Balay } 590*49b5e25fSSatish Balay } 591*49b5e25fSSatish Balay } 592*49b5e25fSSatish Balay PetscFunctionReturn(0); 593*49b5e25fSSatish Balay } 594*49b5e25fSSatish Balay 595*49b5e25fSSatish Balay #undef __FUNC__ 596*49b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw" 597*49b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw(Mat A,Viewer viewer) 598*49b5e25fSSatish Balay { 599*49b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 600*49b5e25fSSatish Balay int ierr; 601*49b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 602*49b5e25fSSatish Balay Draw draw; 603*49b5e25fSSatish Balay PetscTruth isnull; 604*49b5e25fSSatish Balay 605*49b5e25fSSatish Balay PetscFunctionBegin; 606*49b5e25fSSatish Balay 607*49b5e25fSSatish Balay ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 608*49b5e25fSSatish Balay ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 609*49b5e25fSSatish Balay 610*49b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 611*49b5e25fSSatish Balay xr = a->m; yr = a->m; h = yr/10.0; w = xr/10.0; 612*49b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 613*49b5e25fSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 614*49b5e25fSSatish Balay ierr = DrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 615*49b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 616*49b5e25fSSatish Balay PetscFunctionReturn(0); 617*49b5e25fSSatish Balay } 618*49b5e25fSSatish Balay 619*49b5e25fSSatish Balay #undef __FUNC__ 620*49b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ" 621*49b5e25fSSatish Balay int MatView_SeqSBAIJ(Mat A,Viewer viewer) 622*49b5e25fSSatish Balay { 623*49b5e25fSSatish Balay int ierr; 624*49b5e25fSSatish Balay PetscTruth issocket,isascii,isbinary,isdraw; 625*49b5e25fSSatish Balay 626*49b5e25fSSatish Balay PetscFunctionBegin; 627*49b5e25fSSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 628*49b5e25fSSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 629*49b5e25fSSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 630*49b5e25fSSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 631*49b5e25fSSatish Balay if (issocket) { 632*49b5e25fSSatish Balay SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported"); 633*49b5e25fSSatish Balay } else if (isascii){ 634*49b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 635*49b5e25fSSatish Balay } else if (isbinary) { 636*49b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr); 637*49b5e25fSSatish Balay } else if (isdraw) { 638*49b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 639*49b5e25fSSatish Balay } else { 640*49b5e25fSSatish Balay SETERRQ1(1,1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 641*49b5e25fSSatish Balay } 642*49b5e25fSSatish Balay PetscFunctionReturn(0); 643*49b5e25fSSatish Balay } 644*49b5e25fSSatish Balay 645*49b5e25fSSatish Balay 646*49b5e25fSSatish Balay #undef __FUNC__ 647*49b5e25fSSatish Balay #define __FUNC__ "MatGetValues_SeqSBAIJ" 648*49b5e25fSSatish Balay int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 649*49b5e25fSSatish Balay { 650*49b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 651*49b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 652*49b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 653*49b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 654*49b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 655*49b5e25fSSatish Balay 656*49b5e25fSSatish Balay PetscFunctionBegin; 657*49b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 658*49b5e25fSSatish Balay row = im[k]; brow = row/bs; 659*49b5e25fSSatish Balay if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 660*49b5e25fSSatish Balay if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 661*49b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 662*49b5e25fSSatish Balay nrow = ailen[brow]; 663*49b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 664*49b5e25fSSatish Balay if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 665*49b5e25fSSatish Balay if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 666*49b5e25fSSatish Balay col = in[l] ; 667*49b5e25fSSatish Balay bcol = col/bs; 668*49b5e25fSSatish Balay cidx = col%bs; 669*49b5e25fSSatish Balay ridx = row%bs; 670*49b5e25fSSatish Balay high = nrow; 671*49b5e25fSSatish Balay low = 0; /* assume unsorted */ 672*49b5e25fSSatish Balay while (high-low > 5) { 673*49b5e25fSSatish Balay t = (low+high)/2; 674*49b5e25fSSatish Balay if (rp[t] > bcol) high = t; 675*49b5e25fSSatish Balay else low = t; 676*49b5e25fSSatish Balay } 677*49b5e25fSSatish Balay for (i=low; i<high; i++) { 678*49b5e25fSSatish Balay if (rp[i] > bcol) break; 679*49b5e25fSSatish Balay if (rp[i] == bcol) { 680*49b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 681*49b5e25fSSatish Balay goto finished; 682*49b5e25fSSatish Balay } 683*49b5e25fSSatish Balay } 684*49b5e25fSSatish Balay *v++ = zero; 685*49b5e25fSSatish Balay finished:; 686*49b5e25fSSatish Balay } 687*49b5e25fSSatish Balay } 688*49b5e25fSSatish Balay PetscFunctionReturn(0); 689*49b5e25fSSatish Balay } 690*49b5e25fSSatish Balay 691*49b5e25fSSatish Balay 692*49b5e25fSSatish Balay #undef __FUNC__ 693*49b5e25fSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqSBAIJ" 694*49b5e25fSSatish Balay int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 695*49b5e25fSSatish Balay { 696*49b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 697*49b5e25fSSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 698*49b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 699*49b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 700*49b5e25fSSatish Balay Scalar *value = v; 701*49b5e25fSSatish Balay MatScalar *ap,*aa=a->a,*bap; 702*49b5e25fSSatish Balay 703*49b5e25fSSatish Balay PetscFunctionBegin; 704*49b5e25fSSatish Balay if (roworiented) { 705*49b5e25fSSatish Balay stepval = (n-1)*bs; 706*49b5e25fSSatish Balay } else { 707*49b5e25fSSatish Balay stepval = (m-1)*bs; 708*49b5e25fSSatish Balay } 709*49b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 710*49b5e25fSSatish Balay row = im[k]; 711*49b5e25fSSatish Balay if (row < 0) continue; 712*49b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 713*49b5e25fSSatish Balay if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 714*49b5e25fSSatish Balay #endif 715*49b5e25fSSatish Balay rp = aj + ai[row]; 716*49b5e25fSSatish Balay ap = aa + bs2*ai[row]; 717*49b5e25fSSatish Balay rmax = imax[row]; 718*49b5e25fSSatish Balay nrow = ailen[row]; 719*49b5e25fSSatish Balay low = 0; 720*49b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 721*49b5e25fSSatish Balay if (in[l] < 0) continue; 722*49b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 723*49b5e25fSSatish Balay if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 724*49b5e25fSSatish Balay #endif 725*49b5e25fSSatish Balay col = in[l]; 726*49b5e25fSSatish Balay if (roworiented) { 727*49b5e25fSSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 728*49b5e25fSSatish Balay } else { 729*49b5e25fSSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 730*49b5e25fSSatish Balay } 731*49b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 732*49b5e25fSSatish Balay while (high-low > 7) { 733*49b5e25fSSatish Balay t = (low+high)/2; 734*49b5e25fSSatish Balay if (rp[t] > col) high = t; 735*49b5e25fSSatish Balay else low = t; 736*49b5e25fSSatish Balay } 737*49b5e25fSSatish Balay for (i=low; i<high; i++) { 738*49b5e25fSSatish Balay if (rp[i] > col) break; 739*49b5e25fSSatish Balay if (rp[i] == col) { 740*49b5e25fSSatish Balay bap = ap + bs2*i; 741*49b5e25fSSatish Balay if (roworiented) { 742*49b5e25fSSatish Balay if (is == ADD_VALUES) { 743*49b5e25fSSatish Balay for (ii=0; ii<bs; ii++,value+=stepval) { 744*49b5e25fSSatish Balay for (jj=ii; jj<bs2; jj+=bs) { 745*49b5e25fSSatish Balay bap[jj] += *value++; 746*49b5e25fSSatish Balay } 747*49b5e25fSSatish Balay } 748*49b5e25fSSatish Balay } else { 749*49b5e25fSSatish Balay for (ii=0; ii<bs; ii++,value+=stepval) { 750*49b5e25fSSatish Balay for (jj=ii; jj<bs2; jj+=bs) { 751*49b5e25fSSatish Balay bap[jj] = *value++; 752*49b5e25fSSatish Balay } 753*49b5e25fSSatish Balay } 754*49b5e25fSSatish Balay } 755*49b5e25fSSatish Balay } else { 756*49b5e25fSSatish Balay if (is == ADD_VALUES) { 757*49b5e25fSSatish Balay for (ii=0; ii<bs; ii++,value+=stepval) { 758*49b5e25fSSatish Balay for (jj=0; jj<bs; jj++) { 759*49b5e25fSSatish Balay *bap++ += *value++; 760*49b5e25fSSatish Balay } 761*49b5e25fSSatish Balay } 762*49b5e25fSSatish Balay } else { 763*49b5e25fSSatish Balay for (ii=0; ii<bs; ii++,value+=stepval) { 764*49b5e25fSSatish Balay for (jj=0; jj<bs; jj++) { 765*49b5e25fSSatish Balay *bap++ = *value++; 766*49b5e25fSSatish Balay } 767*49b5e25fSSatish Balay } 768*49b5e25fSSatish Balay } 769*49b5e25fSSatish Balay } 770*49b5e25fSSatish Balay goto noinsert2; 771*49b5e25fSSatish Balay } 772*49b5e25fSSatish Balay } 773*49b5e25fSSatish Balay if (nonew == 1) goto noinsert2; 774*49b5e25fSSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 775*49b5e25fSSatish Balay if (nrow >= rmax) { 776*49b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 777*49b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 778*49b5e25fSSatish Balay MatScalar *new_a; 779*49b5e25fSSatish Balay 780*49b5e25fSSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 781*49b5e25fSSatish Balay 782*49b5e25fSSatish Balay /* malloc new storage space */ 783*49b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 784*49b5e25fSSatish Balay new_a = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a); 785*49b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 786*49b5e25fSSatish Balay new_i = new_j + new_nz; 787*49b5e25fSSatish Balay 788*49b5e25fSSatish Balay /* copy over old data into new slots */ 789*49b5e25fSSatish Balay for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 790*49b5e25fSSatish Balay for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 791*49b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 792*49b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[row] - nrow); 793*49b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 794*49b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 795*49b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 796*49b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 797*49b5e25fSSatish Balay /* free up old matrix storage */ 798*49b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 799*49b5e25fSSatish Balay if (!a->singlemalloc) { 800*49b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 801*49b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 802*49b5e25fSSatish Balay } 803*49b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 804*49b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 805*49b5e25fSSatish Balay 806*49b5e25fSSatish Balay rp = aj + ai[row]; ap = aa + bs2*ai[row]; 807*49b5e25fSSatish Balay rmax = imax[row] = imax[row] + CHUNKSIZE; 808*49b5e25fSSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 809*49b5e25fSSatish Balay a->maxnz += bs2*CHUNKSIZE; 810*49b5e25fSSatish Balay a->reallocs++; 811*49b5e25fSSatish Balay a->nz++; 812*49b5e25fSSatish Balay } 813*49b5e25fSSatish Balay N = nrow++ - 1; 814*49b5e25fSSatish Balay /* shift up all the later entries in this row */ 815*49b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 816*49b5e25fSSatish Balay rp[ii+1] = rp[ii]; 817*49b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 818*49b5e25fSSatish Balay } 819*49b5e25fSSatish Balay if (N >= i) { 820*49b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 821*49b5e25fSSatish Balay } 822*49b5e25fSSatish Balay rp[i] = col; 823*49b5e25fSSatish Balay bap = ap + bs2*i; 824*49b5e25fSSatish Balay if (roworiented) { 825*49b5e25fSSatish Balay for (ii=0; ii<bs; ii++,value+=stepval) { 826*49b5e25fSSatish Balay for (jj=ii; jj<bs2; jj+=bs) { 827*49b5e25fSSatish Balay bap[jj] = *value++; 828*49b5e25fSSatish Balay } 829*49b5e25fSSatish Balay } 830*49b5e25fSSatish Balay } else { 831*49b5e25fSSatish Balay for (ii=0; ii<bs; ii++,value+=stepval) { 832*49b5e25fSSatish Balay for (jj=0; jj<bs; jj++) { 833*49b5e25fSSatish Balay *bap++ = *value++; 834*49b5e25fSSatish Balay } 835*49b5e25fSSatish Balay } 836*49b5e25fSSatish Balay } 837*49b5e25fSSatish Balay noinsert2:; 838*49b5e25fSSatish Balay low = i; 839*49b5e25fSSatish Balay } 840*49b5e25fSSatish Balay ailen[row] = nrow; 841*49b5e25fSSatish Balay } 842*49b5e25fSSatish Balay PetscFunctionReturn(0); 843*49b5e25fSSatish Balay } 844*49b5e25fSSatish Balay 845*49b5e25fSSatish Balay #undef __FUNC__ 846*49b5e25fSSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqSBAIJ" 847*49b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 848*49b5e25fSSatish Balay { 849*49b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 850*49b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 851*49b5e25fSSatish Balay int m = a->m,*ip,N,*ailen = a->ilen; 852*49b5e25fSSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 853*49b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 854*49b5e25fSSatish Balay 855*49b5e25fSSatish Balay PetscFunctionBegin; 856*49b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 857*49b5e25fSSatish Balay 858*49b5e25fSSatish Balay if (m) rmax = ailen[0]; 859*49b5e25fSSatish Balay for (i=1; i<mbs; i++) { 860*49b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 861*49b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 862*49b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 863*49b5e25fSSatish Balay if (fshift) { 864*49b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 865*49b5e25fSSatish Balay N = ailen[i]; 866*49b5e25fSSatish Balay for (j=0; j<N; j++) { 867*49b5e25fSSatish Balay ip[j-fshift] = ip[j]; 868*49b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 869*49b5e25fSSatish Balay } 870*49b5e25fSSatish Balay } 871*49b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 872*49b5e25fSSatish Balay } 873*49b5e25fSSatish Balay if (mbs) { 874*49b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 875*49b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 876*49b5e25fSSatish Balay } 877*49b5e25fSSatish Balay /* reset ilen and imax for each row */ 878*49b5e25fSSatish Balay for (i=0; i<mbs; i++) { 879*49b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 880*49b5e25fSSatish Balay } 881*49b5e25fSSatish Balay a->s_nz = ai[mbs]; 882*49b5e25fSSatish Balay 883*49b5e25fSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 884*49b5e25fSSatish Balay if (fshift && a->diag) { 885*49b5e25fSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 886*49b5e25fSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 887*49b5e25fSSatish Balay a->diag = 0; 888*49b5e25fSSatish Balay } 889*49b5e25fSSatish Balay PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 890*49b5e25fSSatish Balay m,a->m,a->bs,fshift*bs2,a->s_nz*bs2); 891*49b5e25fSSatish Balay PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n", 892*49b5e25fSSatish Balay a->reallocs); 893*49b5e25fSSatish Balay PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 894*49b5e25fSSatish Balay a->reallocs = 0; 895*49b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 896*49b5e25fSSatish Balay 897*49b5e25fSSatish Balay PetscFunctionReturn(0); 898*49b5e25fSSatish Balay } 899*49b5e25fSSatish Balay 900*49b5e25fSSatish Balay /* 901*49b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 902*49b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 903*49b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 904*49b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 905*49b5e25fSSatish Balay */ 906*49b5e25fSSatish Balay #undef __FUNC__ 907*49b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 908*49b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 909*49b5e25fSSatish Balay { 910*49b5e25fSSatish Balay int i,j,k,row; 911*49b5e25fSSatish Balay PetscTruth flg; 912*49b5e25fSSatish Balay 913*49b5e25fSSatish Balay PetscFunctionBegin; 914*49b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 915*49b5e25fSSatish Balay row = idx[i]; 916*49b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 917*49b5e25fSSatish Balay sizes[j] = 1; 918*49b5e25fSSatish Balay i++; 919*49b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 920*49b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 921*49b5e25fSSatish Balay i++; 922*49b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 923*49b5e25fSSatish Balay flg = PETSC_TRUE; 924*49b5e25fSSatish Balay for (k=1; k<bs; k++) { 925*49b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 926*49b5e25fSSatish Balay flg = PETSC_FALSE; 927*49b5e25fSSatish Balay break; 928*49b5e25fSSatish Balay } 929*49b5e25fSSatish Balay } 930*49b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 931*49b5e25fSSatish Balay sizes[j] = bs; 932*49b5e25fSSatish Balay i+= bs; 933*49b5e25fSSatish Balay } else { 934*49b5e25fSSatish Balay sizes[j] = 1; 935*49b5e25fSSatish Balay i++; 936*49b5e25fSSatish Balay } 937*49b5e25fSSatish Balay } 938*49b5e25fSSatish Balay } 939*49b5e25fSSatish Balay *bs_max = j; 940*49b5e25fSSatish Balay PetscFunctionReturn(0); 941*49b5e25fSSatish Balay } 942*49b5e25fSSatish Balay 943*49b5e25fSSatish Balay #undef __FUNC__ 944*49b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ" 945*49b5e25fSSatish Balay int MatZeroRows_SeqSBAIJ(Mat A,IS is,Scalar *diag) 946*49b5e25fSSatish Balay { 947*49b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data; 948*49b5e25fSSatish Balay int ierr,i,j,k,count,is_n,*is_idx,*rows; 949*49b5e25fSSatish Balay int bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max; 950*49b5e25fSSatish Balay Scalar zero = 0.0; 951*49b5e25fSSatish Balay MatScalar *aa; 952*49b5e25fSSatish Balay 953*49b5e25fSSatish Balay PetscFunctionBegin; 954*49b5e25fSSatish Balay /* Make a copy of the IS and sort it */ 955*49b5e25fSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 956*49b5e25fSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 957*49b5e25fSSatish Balay 958*49b5e25fSSatish Balay /* allocate memory for rows,sizes */ 959*49b5e25fSSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 960*49b5e25fSSatish Balay sizes = rows + is_n; 961*49b5e25fSSatish Balay 962*49b5e25fSSatish Balay /* initialize copy IS values to rows, and sort them */ 963*49b5e25fSSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 964*49b5e25fSSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 965*49b5e25fSSatish Balay if (sbaij->keepzeroedrows) { /* do not change nonzero structure */ 966*49b5e25fSSatish Balay for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */ 967*49b5e25fSSatish Balay bs_max = is_n; /* bs_max: num. of contiguous block row in the row */ 968*49b5e25fSSatish Balay } else { 969*49b5e25fSSatish Balay ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 970*49b5e25fSSatish Balay } 971*49b5e25fSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 972*49b5e25fSSatish Balay 973*49b5e25fSSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 974*49b5e25fSSatish Balay row = rows[j]; /* row to be zeroed */ 975*49b5e25fSSatish Balay if (row < 0 || row > sbaij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row); 976*49b5e25fSSatish Balay count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */ 977*49b5e25fSSatish Balay aa = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs); 978*49b5e25fSSatish Balay if (sizes[i] == bs && !sbaij->keepzeroedrows) { 979*49b5e25fSSatish Balay if (diag) { 980*49b5e25fSSatish Balay if (sbaij->ilen[row/bs] > 0) { 981*49b5e25fSSatish Balay sbaij->ilen[row/bs] = 1; 982*49b5e25fSSatish Balay sbaij->j[sbaij->i[row/bs]] = row/bs; 983*49b5e25fSSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 984*49b5e25fSSatish Balay } 985*49b5e25fSSatish Balay /* Now insert all the diagoanl values for this bs */ 986*49b5e25fSSatish Balay for (k=0; k<bs; k++) { 987*49b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 988*49b5e25fSSatish Balay } 989*49b5e25fSSatish Balay } else { /* (!diag) */ 990*49b5e25fSSatish Balay sbaij->ilen[row/bs] = 0; 991*49b5e25fSSatish Balay } /* end (!diag) */ 992*49b5e25fSSatish Balay } else { /* (sizes[i] != bs), broken block */ 993*49b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG) 994*49b5e25fSSatish Balay if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1"); 995*49b5e25fSSatish Balay #endif 996*49b5e25fSSatish Balay for (k=0; k<count; k++) { 997*49b5e25fSSatish Balay aa[0] = zero; 998*49b5e25fSSatish Balay aa+=bs; 999*49b5e25fSSatish Balay } 1000*49b5e25fSSatish Balay if (diag) { 1001*49b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 1002*49b5e25fSSatish Balay } 1003*49b5e25fSSatish Balay } 1004*49b5e25fSSatish Balay } 1005*49b5e25fSSatish Balay 1006*49b5e25fSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 1007*49b5e25fSSatish Balay ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1008*49b5e25fSSatish Balay ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1009*49b5e25fSSatish Balay PetscFunctionReturn(0); 1010*49b5e25fSSatish Balay } 1011*49b5e25fSSatish Balay 1012*49b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 1013*49b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 1014*49b5e25fSSatish Balay */ 1015*49b5e25fSSatish Balay 1016*49b5e25fSSatish Balay #undef __FUNC__ 1017*49b5e25fSSatish Balay #define __FUNC__ "MatSetValues_SeqSBAIJ" 1018*49b5e25fSSatish Balay int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 1019*49b5e25fSSatish Balay { 1020*49b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1021*49b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 1022*49b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 1023*49b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 1024*49b5e25fSSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 1025*49b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 1026*49b5e25fSSatish Balay int ambs=a->mbs; 1027*49b5e25fSSatish Balay 1028*49b5e25fSSatish Balay PetscFunctionBegin; 1029*49b5e25fSSatish Balay 1030*49b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 1031*49b5e25fSSatish Balay row = im[k]; /* row number */ 1032*49b5e25fSSatish Balay brow = row/bs; /* block row number */ 1033*49b5e25fSSatish Balay if (row < 0) continue; 1034*49b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 1035*49b5e25fSSatish Balay if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m); 1036*49b5e25fSSatish Balay #endif 1037*49b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 1038*49b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 1039*49b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 1040*49b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 1041*49b5e25fSSatish Balay low = 0; 1042*49b5e25fSSatish Balay 1043*49b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 1044*49b5e25fSSatish Balay if (in[l] < 0) continue; 1045*49b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 1046*49b5e25fSSatish Balay if (in[l] >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->m); 1047*49b5e25fSSatish Balay #endif 1048*49b5e25fSSatish Balay col = in[l]; 1049*49b5e25fSSatish Balay bcol = col/bs; /* block col number */ 1050*49b5e25fSSatish Balay 1051*49b5e25fSSatish Balay if (brow <= bcol){ 1052*49b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 1053*49b5e25fSSatish Balay /* if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ */ 1054*49b5e25fSSatish Balay /* element value a(k,l) */ 1055*49b5e25fSSatish Balay if (roworiented) { 1056*49b5e25fSSatish Balay value = v[l + k*n]; 1057*49b5e25fSSatish Balay } else { 1058*49b5e25fSSatish Balay value = v[k + l*m]; 1059*49b5e25fSSatish Balay } 1060*49b5e25fSSatish Balay 1061*49b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 1062*49b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 1063*49b5e25fSSatish Balay while (high-low > 7) { 1064*49b5e25fSSatish Balay t = (low+high)/2; 1065*49b5e25fSSatish Balay if (rp[t] > bcol) high = t; 1066*49b5e25fSSatish Balay else low = t; 1067*49b5e25fSSatish Balay } 1068*49b5e25fSSatish Balay for (i=low; i<high; i++) { 1069*49b5e25fSSatish Balay /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */ 1070*49b5e25fSSatish Balay if (rp[i] > bcol) break; 1071*49b5e25fSSatish Balay if (rp[i] == bcol) { 1072*49b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 1073*49b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 1074*49b5e25fSSatish Balay else *bap = value; 1075*49b5e25fSSatish Balay goto noinsert1; 1076*49b5e25fSSatish Balay } 1077*49b5e25fSSatish Balay } 1078*49b5e25fSSatish Balay 1079*49b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 1080*49b5e25fSSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 1081*49b5e25fSSatish Balay if (nrow >= rmax) { 1082*49b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 1083*49b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 1084*49b5e25fSSatish Balay MatScalar *new_a; 1085*49b5e25fSSatish Balay 1086*49b5e25fSSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 1087*49b5e25fSSatish Balay 1088*49b5e25fSSatish Balay /* Malloc new storage space */ 1089*49b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 1090*49b5e25fSSatish Balay new_a = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a); 1091*49b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 1092*49b5e25fSSatish Balay new_i = new_j + new_nz; 1093*49b5e25fSSatish Balay 1094*49b5e25fSSatish Balay /* copy over old data into new slots */ 1095*49b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 1096*49b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1097*49b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 1098*49b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1099*49b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1100*49b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1101*49b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1102*49b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 1103*49b5e25fSSatish Balay /* free up old matrix storage */ 1104*49b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1105*49b5e25fSSatish Balay if (!a->singlemalloc) { 1106*49b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1107*49b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1108*49b5e25fSSatish Balay } 1109*49b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 1110*49b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 1111*49b5e25fSSatish Balay 1112*49b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 1113*49b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 1114*49b5e25fSSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 1115*49b5e25fSSatish Balay a->s_maxnz += bs2*CHUNKSIZE; 1116*49b5e25fSSatish Balay a->reallocs++; 1117*49b5e25fSSatish Balay a->s_nz++; 1118*49b5e25fSSatish Balay } 1119*49b5e25fSSatish Balay 1120*49b5e25fSSatish Balay N = nrow++ - 1; 1121*49b5e25fSSatish Balay /* shift up all the later entries in this row */ 1122*49b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 1123*49b5e25fSSatish Balay rp[ii+1] = rp[ii]; 1124*49b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 1125*49b5e25fSSatish Balay } 1126*49b5e25fSSatish Balay if (N>=i) { 1127*49b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1128*49b5e25fSSatish Balay } 1129*49b5e25fSSatish Balay rp[i] = bcol; 1130*49b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 1131*49b5e25fSSatish Balay noinsert1:; 1132*49b5e25fSSatish Balay low = i; 1133*49b5e25fSSatish Balay /* } */ 1134*49b5e25fSSatish Balay } /* end of if .. if.. */ 1135*49b5e25fSSatish Balay } /* end of loop over added columns */ 1136*49b5e25fSSatish Balay ailen[brow] = nrow; 1137*49b5e25fSSatish Balay } /* end of loop over added rows */ 1138*49b5e25fSSatish Balay 1139*49b5e25fSSatish Balay PetscFunctionReturn(0); 1140*49b5e25fSSatish Balay } 1141*49b5e25fSSatish Balay 1142*49b5e25fSSatish Balay extern int MatLUFactorSymbolic_SeqSBAIJ(Mat,IS,IS,MatLUInfo*,Mat*); 1143*49b5e25fSSatish Balay extern int MatLUFactor_SeqSBAIJ(Mat,IS,IS,MatLUInfo*); 1144*49b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int); 1145*49b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 1146*49b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 1147*49b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec); 1148*49b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec); 1149*49b5e25fSSatish Balay extern int MatScale_SeqSBAIJ(Scalar*,Mat); 1150*49b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 1151*49b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 1152*49b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec); 1153*49b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 1154*49b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 1155*49b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat); 1156*49b5e25fSSatish Balay 1157*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 1158*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 1159*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 1160*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 1161*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 1162*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 1163*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 1164*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 1165*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec); 1166*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec); 1167*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec); 1168*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec); 1169*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec); 1170*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec); 1171*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec); 1172*49b5e25fSSatish Balay 1173*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 1174*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 1175*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 1176*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 1177*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 1178*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 1179*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 1180*49b5e25fSSatish Balay 1181*49b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 1182*49b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 1183*49b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 1184*49b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 1185*49b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 1186*49b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 1187*49b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 1188*49b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 1189*49b5e25fSSatish Balay 1190*49b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 1191*49b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 1192*49b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 1193*49b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 1194*49b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 1195*49b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 1196*49b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 1197*49b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 1198*49b5e25fSSatish Balay 1199*49b5e25fSSatish Balay #undef __FUNC__ 1200*49b5e25fSSatish Balay #define __FUNC__ "MatILUFactor_SeqSBAIJ" 1201*49b5e25fSSatish Balay int MatILUFactor_SeqSBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 1202*49b5e25fSSatish Balay { 1203*49b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 1204*49b5e25fSSatish Balay Mat outA; 1205*49b5e25fSSatish Balay int ierr; 1206*49b5e25fSSatish Balay PetscTruth row_identity,col_identity; 1207*49b5e25fSSatish Balay 1208*49b5e25fSSatish Balay PetscFunctionBegin; 1209*49b5e25fSSatish Balay if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU"); 1210*49b5e25fSSatish Balay ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1211*49b5e25fSSatish Balay ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1212*49b5e25fSSatish Balay if (!row_identity || !col_identity) { 1213*49b5e25fSSatish Balay SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU"); 1214*49b5e25fSSatish Balay } 1215*49b5e25fSSatish Balay 1216*49b5e25fSSatish Balay outA = inA; 1217*49b5e25fSSatish Balay inA->factor = FACTOR_LU; 1218*49b5e25fSSatish Balay 1219*49b5e25fSSatish Balay if (!a->diag) { 1220*49b5e25fSSatish Balay ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr); 1221*49b5e25fSSatish Balay } 1222*49b5e25fSSatish Balay /* 1223*49b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 1224*49b5e25fSSatish Balay for ILU(0) factorization with natural ordering 1225*49b5e25fSSatish Balay */ 1226*49b5e25fSSatish Balay switch (a->bs) { 1227*49b5e25fSSatish Balay case 1: 1228*49b5e25fSSatish Balay inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 1229*49b5e25fSSatish Balay PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 1230*49b5e25fSSatish Balay case 2: 1231*49b5e25fSSatish Balay inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 1232*49b5e25fSSatish Balay inA->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 1233*49b5e25fSSatish Balay inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 1234*49b5e25fSSatish Balay PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 1235*49b5e25fSSatish Balay break; 1236*49b5e25fSSatish Balay case 3: 1237*49b5e25fSSatish Balay inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 1238*49b5e25fSSatish Balay inA->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 1239*49b5e25fSSatish Balay inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering; 1240*49b5e25fSSatish Balay PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 1241*49b5e25fSSatish Balay break; 1242*49b5e25fSSatish Balay case 4: 1243*49b5e25fSSatish Balay inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1244*49b5e25fSSatish Balay inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 1245*49b5e25fSSatish Balay inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering; 1246*49b5e25fSSatish Balay PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 1247*49b5e25fSSatish Balay break; 1248*49b5e25fSSatish Balay case 5: 1249*49b5e25fSSatish Balay inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1250*49b5e25fSSatish Balay inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 1251*49b5e25fSSatish Balay inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering; 1252*49b5e25fSSatish Balay PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 1253*49b5e25fSSatish Balay break; 1254*49b5e25fSSatish Balay case 6: 1255*49b5e25fSSatish Balay inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 1256*49b5e25fSSatish Balay inA->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 1257*49b5e25fSSatish Balay inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering; 1258*49b5e25fSSatish Balay PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 1259*49b5e25fSSatish Balay break; 1260*49b5e25fSSatish Balay case 7: 1261*49b5e25fSSatish Balay inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 1262*49b5e25fSSatish Balay inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering; 1263*49b5e25fSSatish Balay inA->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering; 1264*49b5e25fSSatish Balay PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 1265*49b5e25fSSatish Balay break; 1266*49b5e25fSSatish Balay default: 1267*49b5e25fSSatish Balay a->row = row; 1268*49b5e25fSSatish Balay a->col = col; 1269*49b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1270*49b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 1271*49b5e25fSSatish Balay 1272*49b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1273*49b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1274*49b5e25fSSatish Balay PLogObjectParent(inA,a->icol); 1275*49b5e25fSSatish Balay 1276*49b5e25fSSatish Balay if (!a->solve_work) { 1277*49b5e25fSSatish Balay a->solve_work = (Scalar*)PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 1278*49b5e25fSSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 1279*49b5e25fSSatish Balay } 1280*49b5e25fSSatish Balay } 1281*49b5e25fSSatish Balay 1282*49b5e25fSSatish Balay ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1283*49b5e25fSSatish Balay 1284*49b5e25fSSatish Balay PetscFunctionReturn(0); 1285*49b5e25fSSatish Balay } 1286*49b5e25fSSatish Balay #undef __FUNC__ 1287*49b5e25fSSatish Balay #define __FUNC__ "MatPrintHelp_SeqSBAIJ" 1288*49b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A) 1289*49b5e25fSSatish Balay { 1290*49b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 1291*49b5e25fSSatish Balay MPI_Comm comm = A->comm; 1292*49b5e25fSSatish Balay int ierr; 1293*49b5e25fSSatish Balay 1294*49b5e25fSSatish Balay PetscFunctionBegin; 1295*49b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 1296*49b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1297*49b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 1298*49b5e25fSSatish Balay PetscFunctionReturn(0); 1299*49b5e25fSSatish Balay } 1300*49b5e25fSSatish Balay 1301*49b5e25fSSatish Balay EXTERN_C_BEGIN 1302*49b5e25fSSatish Balay #undef __FUNC__ 1303*49b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1304*49b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 1305*49b5e25fSSatish Balay { 1306*49b5e25fSSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 1307*49b5e25fSSatish Balay int i,nz,n; 1308*49b5e25fSSatish Balay 1309*49b5e25fSSatish Balay PetscFunctionBegin; 1310*49b5e25fSSatish Balay nz = baij->maxnz; 1311*49b5e25fSSatish Balay n = baij->n; 1312*49b5e25fSSatish Balay for (i=0; i<nz; i++) { 1313*49b5e25fSSatish Balay baij->j[i] = indices[i]; 1314*49b5e25fSSatish Balay } 1315*49b5e25fSSatish Balay baij->nz = nz; 1316*49b5e25fSSatish Balay for (i=0; i<n; i++) { 1317*49b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 1318*49b5e25fSSatish Balay } 1319*49b5e25fSSatish Balay 1320*49b5e25fSSatish Balay PetscFunctionReturn(0); 1321*49b5e25fSSatish Balay } 1322*49b5e25fSSatish Balay EXTERN_C_END 1323*49b5e25fSSatish Balay 1324*49b5e25fSSatish Balay #undef __FUNC__ 1325*49b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices" 1326*49b5e25fSSatish Balay /*@ 1327*49b5e25fSSatish Balay MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 1328*49b5e25fSSatish Balay in the matrix. 1329*49b5e25fSSatish Balay 1330*49b5e25fSSatish Balay Input Parameters: 1331*49b5e25fSSatish Balay + mat - the SeqBAIJ matrix 1332*49b5e25fSSatish Balay - indices - the column indices 1333*49b5e25fSSatish Balay 1334*49b5e25fSSatish Balay Level: advanced 1335*49b5e25fSSatish Balay 1336*49b5e25fSSatish Balay Notes: 1337*49b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 1338*49b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 1339*49b5e25fSSatish Balay of the MatSetValues() operation. 1340*49b5e25fSSatish Balay 1341*49b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1342*49b5e25fSSatish Balay MatCreateSeqBAIJ(). 1343*49b5e25fSSatish Balay 1344*49b5e25fSSatish Balay MUST be called before any calls to MatSetValues(); 1345*49b5e25fSSatish Balay 1346*49b5e25fSSatish Balay @*/ 1347*49b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 1348*49b5e25fSSatish Balay { 1349*49b5e25fSSatish Balay int ierr,(*f)(Mat,int *); 1350*49b5e25fSSatish Balay 1351*49b5e25fSSatish Balay PetscFunctionBegin; 1352*49b5e25fSSatish Balay PetscValidHeaderSpecific(mat,MAT_COOKIE); 1353*49b5e25fSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 1354*49b5e25fSSatish Balay if (f) { 1355*49b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 1356*49b5e25fSSatish Balay } else { 1357*49b5e25fSSatish Balay SETERRQ(1,1,"Wrong type of matrix to set column indices"); 1358*49b5e25fSSatish Balay } 1359*49b5e25fSSatish Balay PetscFunctionReturn(0); 1360*49b5e25fSSatish Balay } 1361*49b5e25fSSatish Balay 1362*49b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 1363*49b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 1364*49b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 1365*49b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 1366*49b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 1367*49b5e25fSSatish Balay MatMultAdd_SeqSBAIJ_N, 1368*49b5e25fSSatish Balay MatMultTranspose_SeqSBAIJ, 1369*49b5e25fSSatish Balay MatMultTransposeAdd_SeqSBAIJ, 1370*49b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 1371*49b5e25fSSatish Balay 0, 1372*49b5e25fSSatish Balay 0, 1373*49b5e25fSSatish Balay 0, 1374*49b5e25fSSatish Balay MatLUFactor_SeqSBAIJ, 1375*49b5e25fSSatish Balay 0, 1376*49b5e25fSSatish Balay 0, 1377*49b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 1378*49b5e25fSSatish Balay MatGetInfo_SeqSBAIJ, 1379*49b5e25fSSatish Balay MatEqual_SeqSBAIJ, 1380*49b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 1381*49b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 1382*49b5e25fSSatish Balay MatNorm_SeqSBAIJ, 1383*49b5e25fSSatish Balay 0, 1384*49b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 1385*49b5e25fSSatish Balay 0, 1386*49b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 1387*49b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1388*49b5e25fSSatish Balay MatZeroRows_SeqSBAIJ, 1389*49b5e25fSSatish Balay MatLUFactorSymbolic_SeqSBAIJ, 1390*49b5e25fSSatish Balay MatLUFactorNumeric_SeqSBAIJ_N, 1391*49b5e25fSSatish Balay 0, 1392*49b5e25fSSatish Balay 0, 1393*49b5e25fSSatish Balay MatGetSize_SeqSBAIJ, 1394*49b5e25fSSatish Balay MatGetSize_SeqSBAIJ, 1395*49b5e25fSSatish Balay MatGetOwnershipRange_SeqSBAIJ, 1396*49b5e25fSSatish Balay MatILUFactorSymbolic_SeqSBAIJ, 1397*49b5e25fSSatish Balay 0, 1398*49b5e25fSSatish Balay 0, 1399*49b5e25fSSatish Balay 0, 1400*49b5e25fSSatish Balay MatDuplicate_SeqSBAIJ, 1401*49b5e25fSSatish Balay 0, 1402*49b5e25fSSatish Balay 0, 1403*49b5e25fSSatish Balay MatILUFactor_SeqSBAIJ, 1404*49b5e25fSSatish Balay 0, 1405*49b5e25fSSatish Balay 0, 1406*49b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 1407*49b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 1408*49b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 1409*49b5e25fSSatish Balay 0, 1410*49b5e25fSSatish Balay MatPrintHelp_SeqSBAIJ, 1411*49b5e25fSSatish Balay MatScale_SeqSBAIJ, 1412*49b5e25fSSatish Balay 0, 1413*49b5e25fSSatish Balay 0, 1414*49b5e25fSSatish Balay 0, 1415*49b5e25fSSatish Balay MatGetBlockSize_SeqSBAIJ, 1416*49b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 1417*49b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 1418*49b5e25fSSatish Balay 0, 1419*49b5e25fSSatish Balay 0, 1420*49b5e25fSSatish Balay 0, 1421*49b5e25fSSatish Balay 0, 1422*49b5e25fSSatish Balay 0, 1423*49b5e25fSSatish Balay 0, 1424*49b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 1425*49b5e25fSSatish Balay MatGetSubMatrix_SeqSBAIJ, 1426*49b5e25fSSatish Balay 0, 1427*49b5e25fSSatish Balay 0, 1428*49b5e25fSSatish Balay MatGetMaps_Petsc}; 1429*49b5e25fSSatish Balay 1430*49b5e25fSSatish Balay EXTERN_C_BEGIN 1431*49b5e25fSSatish Balay #undef __FUNC__ 1432*49b5e25fSSatish Balay #define __FUNC__ "MatStoreValues_SeqSBAIJ" 1433*49b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat) 1434*49b5e25fSSatish Balay { 1435*49b5e25fSSatish Balay Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1436*49b5e25fSSatish Balay int nz = aij->i[aij->m]*aij->bs*aij->bs2; 1437*49b5e25fSSatish Balay int ierr; 1438*49b5e25fSSatish Balay 1439*49b5e25fSSatish Balay PetscFunctionBegin; 1440*49b5e25fSSatish Balay if (aij->nonew != 1) { 1441*49b5e25fSSatish Balay SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1442*49b5e25fSSatish Balay } 1443*49b5e25fSSatish Balay 1444*49b5e25fSSatish Balay /* allocate space for values if not already there */ 1445*49b5e25fSSatish Balay if (!aij->saved_values) { 1446*49b5e25fSSatish Balay aij->saved_values = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 1447*49b5e25fSSatish Balay } 1448*49b5e25fSSatish Balay 1449*49b5e25fSSatish Balay /* copy values over */ 1450*49b5e25fSSatish Balay ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 1451*49b5e25fSSatish Balay PetscFunctionReturn(0); 1452*49b5e25fSSatish Balay } 1453*49b5e25fSSatish Balay EXTERN_C_END 1454*49b5e25fSSatish Balay 1455*49b5e25fSSatish Balay EXTERN_C_BEGIN 1456*49b5e25fSSatish Balay #undef __FUNC__ 1457*49b5e25fSSatish Balay #define __FUNC__ "MatRetrieveValues_SeqSBAIJ" 1458*49b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat) 1459*49b5e25fSSatish Balay { 1460*49b5e25fSSatish Balay Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1461*49b5e25fSSatish Balay int nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr; 1462*49b5e25fSSatish Balay 1463*49b5e25fSSatish Balay PetscFunctionBegin; 1464*49b5e25fSSatish Balay if (aij->nonew != 1) { 1465*49b5e25fSSatish Balay SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1466*49b5e25fSSatish Balay } 1467*49b5e25fSSatish Balay if (!aij->saved_values) { 1468*49b5e25fSSatish Balay SETERRQ(1,1,"Must call MatStoreValues(A);first"); 1469*49b5e25fSSatish Balay } 1470*49b5e25fSSatish Balay 1471*49b5e25fSSatish Balay /* copy values over */ 1472*49b5e25fSSatish Balay ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 1473*49b5e25fSSatish Balay PetscFunctionReturn(0); 1474*49b5e25fSSatish Balay } 1475*49b5e25fSSatish Balay EXTERN_C_END 1476*49b5e25fSSatish Balay 1477*49b5e25fSSatish Balay /*------------------------------------ New 1 -------------------------------*/ 1478*49b5e25fSSatish Balay #undef __FUNC__ 1479*49b5e25fSSatish Balay #define __FUNC__ "MatCreateSeqSBAIJ" 1480*49b5e25fSSatish Balay /*@C 1481*49b5e25fSSatish Balay MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1482*49b5e25fSSatish Balay compressed row) format. For good matrix assembly performance the 1483*49b5e25fSSatish Balay user should preallocate the matrix storage by setting the parameter nz 1484*49b5e25fSSatish Balay (or the array nnz). By setting these parameters accurately, performance 1485*49b5e25fSSatish Balay during matrix assembly can be increased by more than a factor of 50. 1486*49b5e25fSSatish Balay 1487*49b5e25fSSatish Balay Collective on MPI_Comm 1488*49b5e25fSSatish Balay 1489*49b5e25fSSatish Balay Input Parameters: 1490*49b5e25fSSatish Balay + comm - MPI communicator, set to PETSC_COMM_SELF 1491*49b5e25fSSatish Balay . bs - size of block 1492*49b5e25fSSatish Balay . m - number of rows, or number of columns 1493*49b5e25fSSatish Balay . nz - number of block nonzeros per block row (same for all rows) 1494*49b5e25fSSatish Balay - nnz - array containing the number of block nonzeros in the various block rows 1495*49b5e25fSSatish Balay (possibly different for each block row) or PETSC_NULL 1496*49b5e25fSSatish Balay 1497*49b5e25fSSatish Balay Output Parameter: 1498*49b5e25fSSatish Balay . A - the symmetric matrix 1499*49b5e25fSSatish Balay 1500*49b5e25fSSatish Balay Options Database Keys: 1501*49b5e25fSSatish Balay . -mat_no_unroll - uses code that does not unroll the loops in the 1502*49b5e25fSSatish Balay block calculations (much slower) 1503*49b5e25fSSatish Balay . -mat_block_size - size of the blocks to use 1504*49b5e25fSSatish Balay 1505*49b5e25fSSatish Balay Level: intermediate 1506*49b5e25fSSatish Balay 1507*49b5e25fSSatish Balay Notes: 1508*49b5e25fSSatish Balay The block AIJ format is fully compatible with standard Fortran 77 1509*49b5e25fSSatish Balay storage. That is, the stored row and column indices can begin at 1510*49b5e25fSSatish Balay either one (as in Fortran) or zero. See the users' manual for details. 1511*49b5e25fSSatish Balay 1512*49b5e25fSSatish Balay Specify the preallocated storage with either nz or nnz (not both). 1513*49b5e25fSSatish Balay Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1514*49b5e25fSSatish Balay allocation. For additional details, see the users manual chapter on 1515*49b5e25fSSatish Balay matrices. 1516*49b5e25fSSatish Balay 1517*49b5e25fSSatish Balay .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 1518*49b5e25fSSatish Balay @*/ 1519*49b5e25fSSatish Balay int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1520*49b5e25fSSatish Balay { 1521*49b5e25fSSatish Balay Mat B; 1522*49b5e25fSSatish Balay Mat_SeqSBAIJ *b; 1523*49b5e25fSSatish Balay int i,len,ierr,mbs,nbs,bs2,size; 1524*49b5e25fSSatish Balay PetscTruth flg; 1525*49b5e25fSSatish Balay int s_nz,s_maxnz,*imax; 1526*49b5e25fSSatish Balay 1527*49b5e25fSSatish Balay PetscFunctionBegin; 1528*49b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1529*49b5e25fSSatish Balay if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1530*49b5e25fSSatish Balay 1531*49b5e25fSSatish Balay ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1532*49b5e25fSSatish Balay mbs = m/bs; 1533*49b5e25fSSatish Balay /* nbs = n/bs; */ 1534*49b5e25fSSatish Balay bs2 = bs*bs; 1535*49b5e25fSSatish Balay 1536*49b5e25fSSatish Balay if (mbs*bs!=m) { 1537*49b5e25fSSatish Balay SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 1538*49b5e25fSSatish Balay } 1539*49b5e25fSSatish Balay 1540*49b5e25fSSatish Balay if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz); 1541*49b5e25fSSatish Balay if (nnz) { 1542*49b5e25fSSatish Balay for (i=0; i<mbs; i++) { 1543*49b5e25fSSatish Balay if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 1544*49b5e25fSSatish Balay } 1545*49b5e25fSSatish Balay } 1546*49b5e25fSSatish Balay 1547*49b5e25fSSatish Balay *A = 0; 1548*49b5e25fSSatish Balay PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQSBAIJ,"Mat",comm,MatDestroy,MatView); 1549*49b5e25fSSatish Balay PLogObjectCreate(B); 1550*49b5e25fSSatish Balay B->data = (void*)(b = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(b); 1551*49b5e25fSSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1552*49b5e25fSSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1553*49b5e25fSSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 1554*49b5e25fSSatish Balay if (!flg) { 1555*49b5e25fSSatish Balay switch (bs) { 1556*49b5e25fSSatish Balay case 1: 1557*49b5e25fSSatish Balay B->ops->lufactornumeric = MatLUFactorNumeric_SeqSBAIJ_1; 1558*49b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1559*49b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_1; 1560*49b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 1561*49b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1562*49b5e25fSSatish Balay break; 1563*49b5e25fSSatish Balay case 2: 1564*49b5e25fSSatish Balay B->ops->lufactornumeric = MatLUFactorNumeric_SeqSBAIJ_2; 1565*49b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1566*49b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_2; 1567*49b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 1568*49b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1569*49b5e25fSSatish Balay break; 1570*49b5e25fSSatish Balay case 3: 1571*49b5e25fSSatish Balay B->ops->lufactornumeric = MatLUFactorNumeric_SeqSBAIJ_3; 1572*49b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1573*49b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_3; 1574*49b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 1575*49b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1576*49b5e25fSSatish Balay break; 1577*49b5e25fSSatish Balay case 4: 1578*49b5e25fSSatish Balay B->ops->lufactornumeric = MatLUFactorNumeric_SeqSBAIJ_4; 1579*49b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1580*49b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_4; 1581*49b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 1582*49b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1583*49b5e25fSSatish Balay break; 1584*49b5e25fSSatish Balay case 5: 1585*49b5e25fSSatish Balay B->ops->lufactornumeric = MatLUFactorNumeric_SeqSBAIJ_5; 1586*49b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1587*49b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_5; 1588*49b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 1589*49b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1590*49b5e25fSSatish Balay break; 1591*49b5e25fSSatish Balay case 6: 1592*49b5e25fSSatish Balay B->ops->lufactornumeric = MatLUFactorNumeric_SeqSBAIJ_6; 1593*49b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1594*49b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_6; 1595*49b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 1596*49b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1597*49b5e25fSSatish Balay break; 1598*49b5e25fSSatish Balay case 7: 1599*49b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_7; 1600*49b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1601*49b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_7; 1602*49b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1603*49b5e25fSSatish Balay break; 1604*49b5e25fSSatish Balay } 1605*49b5e25fSSatish Balay } 1606*49b5e25fSSatish Balay B->ops->destroy = MatDestroy_SeqSBAIJ; 1607*49b5e25fSSatish Balay B->ops->view = MatView_SeqSBAIJ; 1608*49b5e25fSSatish Balay B->factor = 0; 1609*49b5e25fSSatish Balay B->lupivotthreshold = 1.0; 1610*49b5e25fSSatish Balay B->mapping = 0; 1611*49b5e25fSSatish Balay b->row = 0; 1612*49b5e25fSSatish Balay b->col = 0; 1613*49b5e25fSSatish Balay b->icol = 0; 1614*49b5e25fSSatish Balay b->reallocs = 0; 1615*49b5e25fSSatish Balay b->saved_values = 0; 1616*49b5e25fSSatish Balay 1617*49b5e25fSSatish Balay b->m = m; B->m = m; B->M = m; 1618*49b5e25fSSatish Balay /* b->n = n;*/ 1619*49b5e25fSSatish Balay B->n = m; B->N = m; 1620*49b5e25fSSatish Balay 1621*49b5e25fSSatish Balay ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 1622*49b5e25fSSatish Balay ierr = MapCreateMPI(comm,m,m,&B->cmap);CHKERRQ(ierr); 1623*49b5e25fSSatish Balay 1624*49b5e25fSSatish Balay b->mbs = mbs; 1625*49b5e25fSSatish Balay /* b->nbs = nbs; */ 1626*49b5e25fSSatish Balay b->imax = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(b->imax); 1627*49b5e25fSSatish Balay if (!nnz) { 1628*49b5e25fSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 1629*49b5e25fSSatish Balay else if (nz <= 0) nz = 1; 1630*49b5e25fSSatish Balay for (i=0; i<mbs; i++) { 1631*49b5e25fSSatish Balay b->imax[i] = (nz+1)/2; 1632*49b5e25fSSatish Balay } 1633*49b5e25fSSatish Balay nz = nz*mbs; 1634*49b5e25fSSatish Balay } else { 1635*49b5e25fSSatish Balay nz = 0; 1636*49b5e25fSSatish Balay for (i=0; i<mbs; i++) {b->imax[i] = (nnz[i]+1)/2; nz += nnz[i];} 1637*49b5e25fSSatish Balay } 1638*49b5e25fSSatish Balay s_nz=(nz+mbs)/2; /* total diagonal and superdiagonal nonzero blocks */ 1639*49b5e25fSSatish Balay 1640*49b5e25fSSatish Balay /* allocate the matrix space */ 1641*49b5e25fSSatish Balay len = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int); 1642*49b5e25fSSatish Balay b->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(b->a); 1643*49b5e25fSSatish Balay ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1644*49b5e25fSSatish Balay b->j = (int*)(b->a + s_nz*bs2); 1645*49b5e25fSSatish Balay ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr); 1646*49b5e25fSSatish Balay b->i = b->j + s_nz; 1647*49b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 1648*49b5e25fSSatish Balay 1649*49b5e25fSSatish Balay /* pointer to beginning of each row */ 1650*49b5e25fSSatish Balay b->i[0] = 0; 1651*49b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 1652*49b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 1653*49b5e25fSSatish Balay } 1654*49b5e25fSSatish Balay 1655*49b5e25fSSatish Balay 1656*49b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 1657*49b5e25fSSatish Balay b->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int)); 1658*49b5e25fSSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 1659*49b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 1660*49b5e25fSSatish Balay 1661*49b5e25fSSatish Balay b->bs = bs; 1662*49b5e25fSSatish Balay b->bs2 = bs2; 1663*49b5e25fSSatish Balay b->mbs = mbs; 1664*49b5e25fSSatish Balay b->s_nz = 0; 1665*49b5e25fSSatish Balay b->s_maxnz = s_nz*bs2; 1666*49b5e25fSSatish Balay b->sorted = PETSC_FALSE; 1667*49b5e25fSSatish Balay b->roworiented = PETSC_TRUE; 1668*49b5e25fSSatish Balay b->nonew = 0; 1669*49b5e25fSSatish Balay b->diag = 0; 1670*49b5e25fSSatish Balay b->solve_work = 0; 1671*49b5e25fSSatish Balay b->mult_work = 0; 1672*49b5e25fSSatish Balay b->spptr = 0; 1673*49b5e25fSSatish Balay B->info.nz_unneeded = (PetscReal)b->s_maxnz; 1674*49b5e25fSSatish Balay b->keepzeroedrows = PETSC_FALSE; 1675*49b5e25fSSatish Balay 1676*49b5e25fSSatish Balay *A = B; 1677*49b5e25fSSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 1678*49b5e25fSSatish Balay if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 1679*49b5e25fSSatish Balay 1680*49b5e25fSSatish Balay ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1681*49b5e25fSSatish Balay "MatStoreValues_SeqBAIJ", 1682*49b5e25fSSatish Balay (void*)MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 1683*49b5e25fSSatish Balay ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1684*49b5e25fSSatish Balay "MatRetrieveValues_SeqBAIJ", 1685*49b5e25fSSatish Balay (void*)MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 1686*49b5e25fSSatish Balay ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 1687*49b5e25fSSatish Balay "MatSeqBAIJSetColumnIndices_SeqBAIJ", 1688*49b5e25fSSatish Balay (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 1689*49b5e25fSSatish Balay /* printf("In MatCreateSeqSBAIJ, \n"); 1690*49b5e25fSSatish Balay for (i=0; i<mbs; i++){ 1691*49b5e25fSSatish Balay printf("imax[%d]= %d, ilen[%d]= %d\n", i,b->imax[i], i,b->ilen[i]); 1692*49b5e25fSSatish Balay } */ 1693*49b5e25fSSatish Balay 1694*49b5e25fSSatish Balay PetscFunctionReturn(0); 1695*49b5e25fSSatish Balay } 1696*49b5e25fSSatish Balay 1697*49b5e25fSSatish Balay #undef __FUNC__ 1698*49b5e25fSSatish Balay #define __FUNC__ "MatDuplicate_SeqSBAIJ" 1699*49b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 1700*49b5e25fSSatish Balay { 1701*49b5e25fSSatish Balay Mat C; 1702*49b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 1703*49b5e25fSSatish Balay int i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr; 1704*49b5e25fSSatish Balay 1705*49b5e25fSSatish Balay PetscFunctionBegin; 1706*49b5e25fSSatish Balay if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 1707*49b5e25fSSatish Balay 1708*49b5e25fSSatish Balay *B = 0; 1709*49b5e25fSSatish Balay PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQSBAIJ,"Mat",A->comm,MatDestroy,MatView); 1710*49b5e25fSSatish Balay PLogObjectCreate(C); 1711*49b5e25fSSatish Balay C->data = (void*)(c = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(c); 1712*49b5e25fSSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1713*49b5e25fSSatish Balay C->ops->destroy = MatDestroy_SeqSBAIJ; 1714*49b5e25fSSatish Balay C->ops->view = MatView_SeqSBAIJ; 1715*49b5e25fSSatish Balay C->factor = A->factor; 1716*49b5e25fSSatish Balay c->row = 0; 1717*49b5e25fSSatish Balay c->col = 0; 1718*49b5e25fSSatish Balay c->icol = 0; 1719*49b5e25fSSatish Balay c->saved_values = 0; 1720*49b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 1721*49b5e25fSSatish Balay C->assembled = PETSC_TRUE; 1722*49b5e25fSSatish Balay 1723*49b5e25fSSatish Balay c->m = C->m = a->m; 1724*49b5e25fSSatish Balay c->n = C->n = a->n; 1725*49b5e25fSSatish Balay C->M = a->m; 1726*49b5e25fSSatish Balay C->N = a->n; 1727*49b5e25fSSatish Balay 1728*49b5e25fSSatish Balay c->bs = a->bs; 1729*49b5e25fSSatish Balay c->bs2 = a->bs2; 1730*49b5e25fSSatish Balay c->mbs = a->mbs; 1731*49b5e25fSSatish Balay c->nbs = a->nbs; 1732*49b5e25fSSatish Balay 1733*49b5e25fSSatish Balay c->imax = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 1734*49b5e25fSSatish Balay c->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 1735*49b5e25fSSatish Balay for (i=0; i<mbs; i++) { 1736*49b5e25fSSatish Balay c->imax[i] = a->imax[i]; 1737*49b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 1738*49b5e25fSSatish Balay } 1739*49b5e25fSSatish Balay 1740*49b5e25fSSatish Balay /* allocate the matrix space */ 1741*49b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 1742*49b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 1743*49b5e25fSSatish Balay c->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(c->a); 1744*49b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 1745*49b5e25fSSatish Balay c->i = c->j + nz; 1746*49b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1747*49b5e25fSSatish Balay if (mbs > 0) { 1748*49b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 1749*49b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 1750*49b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 1751*49b5e25fSSatish Balay } else { 1752*49b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 1753*49b5e25fSSatish Balay } 1754*49b5e25fSSatish Balay } 1755*49b5e25fSSatish Balay 1756*49b5e25fSSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 1757*49b5e25fSSatish Balay c->sorted = a->sorted; 1758*49b5e25fSSatish Balay c->roworiented = a->roworiented; 1759*49b5e25fSSatish Balay c->nonew = a->nonew; 1760*49b5e25fSSatish Balay 1761*49b5e25fSSatish Balay if (a->diag) { 1762*49b5e25fSSatish Balay c->diag = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->diag); 1763*49b5e25fSSatish Balay PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1764*49b5e25fSSatish Balay for (i=0; i<mbs; i++) { 1765*49b5e25fSSatish Balay c->diag[i] = a->diag[i]; 1766*49b5e25fSSatish Balay } 1767*49b5e25fSSatish Balay } else c->diag = 0; 1768*49b5e25fSSatish Balay c->s_nz = a->s_nz; 1769*49b5e25fSSatish Balay c->s_maxnz = a->s_maxnz; 1770*49b5e25fSSatish Balay c->solve_work = 0; 1771*49b5e25fSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1772*49b5e25fSSatish Balay c->mult_work = 0; 1773*49b5e25fSSatish Balay *B = C; 1774*49b5e25fSSatish Balay ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 1775*49b5e25fSSatish Balay PetscFunctionReturn(0); 1776*49b5e25fSSatish Balay } 1777*49b5e25fSSatish Balay 1778*49b5e25fSSatish Balay #undef __FUNC__ 1779*49b5e25fSSatish Balay #define __FUNC__ "MatLoad_SeqSBAIJ" 1780*49b5e25fSSatish Balay int MatLoad_SeqSBAIJ(Viewer viewer,MatType type,Mat *A) 1781*49b5e25fSSatish Balay { 1782*49b5e25fSSatish Balay Mat_SeqSBAIJ *a; 1783*49b5e25fSSatish Balay Mat B; 1784*49b5e25fSSatish Balay int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1785*49b5e25fSSatish Balay int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 1786*49b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1787*49b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 1788*49b5e25fSSatish Balay Scalar *aa; 1789*49b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 1790*49b5e25fSSatish Balay 1791*49b5e25fSSatish Balay PetscFunctionBegin; 1792*49b5e25fSSatish Balay ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1793*49b5e25fSSatish Balay bs2 = bs*bs; 1794*49b5e25fSSatish Balay 1795*49b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1796*49b5e25fSSatish Balay if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 1797*49b5e25fSSatish Balay ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 1798*49b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1799*49b5e25fSSatish Balay if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 1800*49b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 1801*49b5e25fSSatish Balay 1802*49b5e25fSSatish Balay if (header[3] < 0) { 1803*49b5e25fSSatish Balay SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqSBAIJ"); 1804*49b5e25fSSatish Balay } 1805*49b5e25fSSatish Balay 1806*49b5e25fSSatish Balay if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 1807*49b5e25fSSatish Balay 1808*49b5e25fSSatish Balay /* 1809*49b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 1810*49b5e25fSSatish Balay divisible by the blocksize 1811*49b5e25fSSatish Balay */ 1812*49b5e25fSSatish Balay mbs = M/bs; 1813*49b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 1814*49b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 1815*49b5e25fSSatish Balay else mbs++; 1816*49b5e25fSSatish Balay if (extra_rows) { 1817*49b5e25fSSatish Balay PLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 1818*49b5e25fSSatish Balay } 1819*49b5e25fSSatish Balay 1820*49b5e25fSSatish Balay /* read in row lengths */ 1821*49b5e25fSSatish Balay rowlengths = (int*)PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 1822*49b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 1823*49b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 1824*49b5e25fSSatish Balay 1825*49b5e25fSSatish Balay /* read in column indices */ 1826*49b5e25fSSatish Balay jj = (int*)PetscMalloc((nz+extra_rows)*sizeof(int));CHKPTRQ(jj); 1827*49b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 1828*49b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 1829*49b5e25fSSatish Balay 1830*49b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 1831*49b5e25fSSatish Balay browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1832*49b5e25fSSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 1833*49b5e25fSSatish Balay mask = (int*)PetscMalloc(2*mbs*sizeof(int));CHKPTRQ(mask); 1834*49b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 1835*49b5e25fSSatish Balay masked = mask + mbs; 1836*49b5e25fSSatish Balay rowcount = 0; nzcount = 0; 1837*49b5e25fSSatish Balay for (i=0; i<mbs; i++) { 1838*49b5e25fSSatish Balay nmask = 0; 1839*49b5e25fSSatish Balay for (j=0; j<bs; j++) { 1840*49b5e25fSSatish Balay kmax = rowlengths[rowcount]; 1841*49b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1842*49b5e25fSSatish Balay tmp = jj[nzcount++]/bs; 1843*49b5e25fSSatish Balay if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1844*49b5e25fSSatish Balay } 1845*49b5e25fSSatish Balay rowcount++; 1846*49b5e25fSSatish Balay } 1847*49b5e25fSSatish Balay browlengths[i] += nmask; 1848*49b5e25fSSatish Balay /* zero out the mask elements we set */ 1849*49b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1850*49b5e25fSSatish Balay } 1851*49b5e25fSSatish Balay 1852*49b5e25fSSatish Balay /* create our matrix */ 1853*49b5e25fSSatish Balay ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 1854*49b5e25fSSatish Balay B = *A; 1855*49b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 1856*49b5e25fSSatish Balay 1857*49b5e25fSSatish Balay /* set matrix "i" values */ 1858*49b5e25fSSatish Balay a->i[0] = 0; 1859*49b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1860*49b5e25fSSatish Balay a->i[i] = a->i[i-1] + browlengths[i-1]; 1861*49b5e25fSSatish Balay a->ilen[i-1] = browlengths[i-1]; 1862*49b5e25fSSatish Balay } 1863*49b5e25fSSatish Balay a->s_nz = 0; 1864*49b5e25fSSatish Balay for (i=0; i<mbs; i++) a->s_nz += browlengths[i]; 1865*49b5e25fSSatish Balay 1866*49b5e25fSSatish Balay /* read in nonzero values */ 1867*49b5e25fSSatish Balay aa = (Scalar*)PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 1868*49b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 1869*49b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 1870*49b5e25fSSatish Balay 1871*49b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 1872*49b5e25fSSatish Balay nzcount = 0; jcount = 0; 1873*49b5e25fSSatish Balay for (i=0; i<mbs; i++) { 1874*49b5e25fSSatish Balay nzcountb = nzcount; 1875*49b5e25fSSatish Balay nmask = 0; 1876*49b5e25fSSatish Balay for (j=0; j<bs; j++) { 1877*49b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 1878*49b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1879*49b5e25fSSatish Balay tmp = jj[nzcount++]/bs; 1880*49b5e25fSSatish Balay if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1881*49b5e25fSSatish Balay } 1882*49b5e25fSSatish Balay rowcount++; 1883*49b5e25fSSatish Balay } 1884*49b5e25fSSatish Balay /* sort the masked values */ 1885*49b5e25fSSatish Balay ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 1886*49b5e25fSSatish Balay 1887*49b5e25fSSatish Balay /* set "j" values into matrix */ 1888*49b5e25fSSatish Balay maskcount = 1; 1889*49b5e25fSSatish Balay for (j=0; j<nmask; j++) { 1890*49b5e25fSSatish Balay a->j[jcount++] = masked[j]; 1891*49b5e25fSSatish Balay mask[masked[j]] = maskcount++; 1892*49b5e25fSSatish Balay } 1893*49b5e25fSSatish Balay /* set "a" values into matrix */ 1894*49b5e25fSSatish Balay ishift = bs2*a->i[i]; 1895*49b5e25fSSatish Balay for (j=0; j<bs; j++) { 1896*49b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 1897*49b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1898*49b5e25fSSatish Balay tmp = jj[nzcountb]/bs ; 1899*49b5e25fSSatish Balay block = mask[tmp] - 1; 1900*49b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 1901*49b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1902*49b5e25fSSatish Balay a->a[idx] = aa[nzcountb++]; 1903*49b5e25fSSatish Balay } 1904*49b5e25fSSatish Balay } 1905*49b5e25fSSatish Balay /* zero out the mask elements we set */ 1906*49b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1907*49b5e25fSSatish Balay } 1908*49b5e25fSSatish Balay if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1909*49b5e25fSSatish Balay 1910*49b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1911*49b5e25fSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1912*49b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 1913*49b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 1914*49b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 1915*49b5e25fSSatish Balay 1916*49b5e25fSSatish Balay B->assembled = PETSC_TRUE; 1917*49b5e25fSSatish Balay 1918*49b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 1919*49b5e25fSSatish Balay PetscFunctionReturn(0); 1920*49b5e25fSSatish Balay } 1921*49b5e25fSSatish Balay 1922*49b5e25fSSatish Balay 1923*49b5e25fSSatish Balay 1924