1421e10b8SBarry Smith #define PETSCMAT_DLL 2421e10b8SBarry Smith 3421e10b8SBarry Smith /* 4421e10b8SBarry Smith This provides a matrix that consists of Mats 5421e10b8SBarry Smith */ 6421e10b8SBarry Smith 7421e10b8SBarry Smith #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 8421e10b8SBarry Smith #include "src/mat/impls/baij/seq/baij.h" /* use the common AIJ data-structure */ 9ccb205f8SBarry Smith #include "petscksp.h" 10421e10b8SBarry Smith 11421e10b8SBarry Smith #define CHUNKSIZE 15 12421e10b8SBarry Smith 13421e10b8SBarry Smith typedef struct { 14421e10b8SBarry Smith SEQAIJHEADER(Mat); 15421e10b8SBarry Smith SEQBAIJHEADER; 16ccb205f8SBarry Smith Mat *diags; 17421e10b8SBarry Smith 186e63c7a1SBarry Smith Vec left,right,middle,workb; /* dummy vectors to perform local parts of product */ 19421e10b8SBarry Smith } Mat_BlockMat; 20421e10b8SBarry Smith 21421e10b8SBarry Smith #undef __FUNCT__ 22290bbb0aSBarry Smith #define __FUNCT__ "MatRelax_BlockMat_Symmetric" 23290bbb0aSBarry Smith PetscErrorCode MatRelax_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 24290bbb0aSBarry Smith { 25290bbb0aSBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 26290bbb0aSBarry Smith PetscScalar *x; 27290bbb0aSBarry Smith const Mat *v = a->a; 28290bbb0aSBarry Smith const PetscScalar *b; 29290bbb0aSBarry Smith PetscErrorCode ierr; 30290bbb0aSBarry Smith PetscInt n = A->cmap.n,i,mbs = n/A->rmap.bs,j,bs = A->rmap.bs; 31290bbb0aSBarry Smith const PetscInt *idx; 32290bbb0aSBarry Smith IS row,col; 33290bbb0aSBarry Smith MatFactorInfo info; 346e63c7a1SBarry Smith Vec left = a->left,right = a->right, middle = a->middle; 35290bbb0aSBarry Smith Mat *diag; 36290bbb0aSBarry Smith 37290bbb0aSBarry Smith PetscFunctionBegin; 38290bbb0aSBarry Smith its = its*lits; 39290bbb0aSBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 40290bbb0aSBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat"); 41290bbb0aSBarry Smith if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 42290bbb0aSBarry Smith if (fshift) SETERRQ(PETSC_ERR_SUP,"No support yet for fshift"); 43290bbb0aSBarry Smith if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) 44290bbb0aSBarry Smith SETERRQ(PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep"); 45290bbb0aSBarry Smith 46290bbb0aSBarry Smith if (!a->diags) { 47290bbb0aSBarry Smith ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr); 48290bbb0aSBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 49290bbb0aSBarry Smith for (i=0; i<mbs; i++) { 50290bbb0aSBarry Smith ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERING_ND,&row,&col);CHKERRQ(ierr); 516e63c7a1SBarry Smith ierr = MatCholeskyFactorSymbolic(a->a[a->diag[i]],row,&info,a->diags+i);CHKERRQ(ierr); 526e63c7a1SBarry Smith ierr = MatCholeskyFactorNumeric(a->a[a->diag[i]],&info,a->diags+i);CHKERRQ(ierr); 53290bbb0aSBarry Smith ierr = ISDestroy(row);CHKERRQ(ierr); 54290bbb0aSBarry Smith ierr = ISDestroy(col);CHKERRQ(ierr); 55290bbb0aSBarry Smith } 566e63c7a1SBarry Smith ierr = VecDuplicate(bb,&a->workb);CHKERRQ(ierr); 57290bbb0aSBarry Smith } 58290bbb0aSBarry Smith diag = a->diags; 59290bbb0aSBarry Smith 60290bbb0aSBarry Smith ierr = VecSet(xx,0.0);CHKERRQ(ierr); 61290bbb0aSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 62290bbb0aSBarry Smith /* copy right hand side because it must be modified during iteration */ 636e63c7a1SBarry Smith ierr = VecCopy(bb,a->workb);CHKERRQ(ierr); 646e63c7a1SBarry Smith ierr = VecGetArray(a->workb,(PetscScalar**)&b);CHKERRQ(ierr); 65290bbb0aSBarry Smith 66290bbb0aSBarry Smith /* need to add code for when initial guess is zero, see MatRelax_SeqAIJ */ 67290bbb0aSBarry Smith while (its--) { 68290bbb0aSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 69290bbb0aSBarry Smith 70290bbb0aSBarry Smith for (i=0; i<mbs; i++) { 716e63c7a1SBarry Smith n = a->i[i+1] - a->i[i] - 1; 726e63c7a1SBarry Smith idx = a->j + a->i[i] + 1; 736e63c7a1SBarry Smith v = a->a + a->i[i] + 1; 74290bbb0aSBarry Smith 75290bbb0aSBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 76290bbb0aSBarry Smith for (j=0; j<n; j++) { 77290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 78290bbb0aSBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 79290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 80290bbb0aSBarry Smith } 81290bbb0aSBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 82290bbb0aSBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 83290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 84290bbb0aSBarry Smith 85290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 86290bbb0aSBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 876e63c7a1SBarry Smith 886e63c7a1SBarry Smith /* now adjust right hand side, see MatRelax_SeqSBAIJ */ 896e63c7a1SBarry Smith for (j=0; j<n; j++) { 906e63c7a1SBarry Smith ierr = MatMultTranspose(v[j],right,left);CHKERRQ(ierr); 916e63c7a1SBarry Smith ierr = VecPlaceArray(middle,b + idx[j]*bs);CHKERRQ(ierr); 926e63c7a1SBarry Smith ierr = VecAXPY(middle,-1.0,left);CHKERRQ(ierr); 936e63c7a1SBarry Smith ierr = VecResetArray(middle);CHKERRQ(ierr); 946e63c7a1SBarry Smith } 95290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 966e63c7a1SBarry Smith 97290bbb0aSBarry Smith } 98290bbb0aSBarry Smith } 99290bbb0aSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 100290bbb0aSBarry Smith 101290bbb0aSBarry Smith for (i=mbs-1; i>=0; i--) { 1026e63c7a1SBarry Smith n = a->i[i+1] - a->i[i] - 1; 1036e63c7a1SBarry Smith idx = a->j + a->i[i] + 1; 1046e63c7a1SBarry Smith v = a->a + a->i[i] + 1; 105290bbb0aSBarry Smith 106290bbb0aSBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 107290bbb0aSBarry Smith for (j=0; j<n; j++) { 108290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 109290bbb0aSBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 110290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 111290bbb0aSBarry Smith } 112290bbb0aSBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 113290bbb0aSBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 114290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 115290bbb0aSBarry Smith 116290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 117290bbb0aSBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 118290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 119290bbb0aSBarry Smith 120290bbb0aSBarry Smith } 121290bbb0aSBarry Smith } 122290bbb0aSBarry Smith } 123290bbb0aSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1246e63c7a1SBarry Smith ierr = VecRestoreArray(a->workb,(PetscScalar**)&b);CHKERRQ(ierr); 125290bbb0aSBarry Smith PetscFunctionReturn(0); 126290bbb0aSBarry Smith } 127290bbb0aSBarry Smith 128290bbb0aSBarry Smith #undef __FUNCT__ 129ccb205f8SBarry Smith #define __FUNCT__ "MatRelax_BlockMat" 130ccb205f8SBarry Smith PetscErrorCode MatRelax_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 131ccb205f8SBarry Smith { 132ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 133ccb205f8SBarry Smith PetscScalar *x; 134ccb205f8SBarry Smith const Mat *v = a->a; 135ccb205f8SBarry Smith const PetscScalar *b; 136ccb205f8SBarry Smith PetscErrorCode ierr; 137ccb205f8SBarry Smith PetscInt n = A->cmap.n,i,mbs = n/A->rmap.bs,j,bs = A->rmap.bs; 138ccb205f8SBarry Smith const PetscInt *idx; 139ccb205f8SBarry Smith IS row,col; 140ccb205f8SBarry Smith MatFactorInfo info; 141ccb205f8SBarry Smith Vec left = a->left,right = a->right; 142ccb205f8SBarry Smith Mat *diag; 143ccb205f8SBarry Smith 144ccb205f8SBarry Smith PetscFunctionBegin; 145ccb205f8SBarry Smith its = its*lits; 146ccb205f8SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 147ccb205f8SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat"); 148ccb205f8SBarry Smith if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 149ccb205f8SBarry Smith if (fshift) SETERRQ(PETSC_ERR_SUP,"No support yet for fshift"); 150ccb205f8SBarry Smith 151ccb205f8SBarry Smith if (!a->diags) { 152ccb205f8SBarry Smith ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr); 153ccb205f8SBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 154ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 155ccb205f8SBarry Smith ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERING_ND,&row,&col);CHKERRQ(ierr); 156ccb205f8SBarry Smith ierr = MatLUFactorSymbolic(a->a[a->diag[i]],row,col,&info,a->diags+i);CHKERRQ(ierr); 157ccb205f8SBarry Smith ierr = MatLUFactorNumeric(a->a[a->diag[i]],&info,a->diags+i);CHKERRQ(ierr); 158ccb205f8SBarry Smith ierr = ISDestroy(row);CHKERRQ(ierr); 159ccb205f8SBarry Smith ierr = ISDestroy(col);CHKERRQ(ierr); 160ccb205f8SBarry Smith } 161ccb205f8SBarry Smith } 162ccb205f8SBarry Smith diag = a->diags; 163ccb205f8SBarry Smith 164ccb205f8SBarry Smith ierr = VecSet(xx,0.0);CHKERRQ(ierr); 165ccb205f8SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 166ccb205f8SBarry Smith ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 167ccb205f8SBarry Smith 168ccb205f8SBarry Smith /* need to add code for when initial guess is zero, see MatRelax_SeqAIJ */ 169ccb205f8SBarry Smith while (its--) { 170ccb205f8SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 171ccb205f8SBarry Smith 172ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 173ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 174ccb205f8SBarry Smith idx = a->j + a->i[i]; 175ccb205f8SBarry Smith v = a->a + a->i[i]; 176ccb205f8SBarry Smith 177ccb205f8SBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 178ccb205f8SBarry Smith for (j=0; j<n; j++) { 179ccb205f8SBarry Smith if (idx[j] != i) { 180ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 181ccb205f8SBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 182ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 183ccb205f8SBarry Smith } 184ccb205f8SBarry Smith } 185ccb205f8SBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 186ccb205f8SBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 187ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 188ccb205f8SBarry Smith 189ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 190ccb205f8SBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 191ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 192ccb205f8SBarry Smith } 193ccb205f8SBarry Smith } 194ccb205f8SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 195ccb205f8SBarry Smith 196ccb205f8SBarry Smith for (i=mbs-1; i>=0; i--) { 197ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 198ccb205f8SBarry Smith idx = a->j + a->i[i]; 199ccb205f8SBarry Smith v = a->a + a->i[i]; 200ccb205f8SBarry Smith 201ccb205f8SBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 202ccb205f8SBarry Smith for (j=0; j<n; j++) { 203ccb205f8SBarry Smith if (idx[j] != i) { 204ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 205ccb205f8SBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 206ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 207ccb205f8SBarry Smith } 208ccb205f8SBarry Smith } 209ccb205f8SBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 210ccb205f8SBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 211ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 212ccb205f8SBarry Smith 213ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 214ccb205f8SBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 215ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 216ccb205f8SBarry Smith 217ccb205f8SBarry Smith } 218ccb205f8SBarry Smith } 219ccb205f8SBarry Smith } 220ccb205f8SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 221ccb205f8SBarry Smith ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 222ccb205f8SBarry Smith PetscFunctionReturn(0); 223ccb205f8SBarry Smith } 224ccb205f8SBarry Smith 225ccb205f8SBarry Smith #undef __FUNCT__ 226421e10b8SBarry Smith #define __FUNCT__ "MatSetValues_BlockMat" 227421e10b8SBarry Smith PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 228421e10b8SBarry Smith { 229421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 230421e10b8SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 231421e10b8SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 232421e10b8SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol; 233421e10b8SBarry Smith PetscErrorCode ierr; 234421e10b8SBarry Smith PetscInt ridx,cidx; 235421e10b8SBarry Smith PetscTruth roworiented=a->roworiented; 236421e10b8SBarry Smith MatScalar value; 237421e10b8SBarry Smith Mat *ap,*aa = a->a; 238421e10b8SBarry Smith 239421e10b8SBarry Smith PetscFunctionBegin; 240421e10b8SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 241421e10b8SBarry Smith row = im[k]; 242421e10b8SBarry Smith brow = row/bs; 243421e10b8SBarry Smith if (row < 0) continue; 244421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 245421e10b8SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 246421e10b8SBarry Smith #endif 247421e10b8SBarry Smith rp = aj + ai[brow]; 248421e10b8SBarry Smith ap = aa + ai[brow]; 249421e10b8SBarry Smith rmax = imax[brow]; 250421e10b8SBarry Smith nrow = ailen[brow]; 251421e10b8SBarry Smith low = 0; 252421e10b8SBarry Smith high = nrow; 253421e10b8SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 254421e10b8SBarry Smith if (in[l] < 0) continue; 255421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 256421e10b8SBarry Smith if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1); 257421e10b8SBarry Smith #endif 258421e10b8SBarry Smith col = in[l]; bcol = col/bs; 2596e63c7a1SBarry Smith if (A->symmetric && brow > bcol) continue; 260421e10b8SBarry Smith ridx = row % bs; cidx = col % bs; 261421e10b8SBarry Smith if (roworiented) { 262421e10b8SBarry Smith value = v[l + k*n]; 263421e10b8SBarry Smith } else { 264421e10b8SBarry Smith value = v[k + l*m]; 265421e10b8SBarry Smith } 266*a0e1a203SBarry Smith /* printf(" brow %d bcol %d\n",brow,bcol);*/ 267421e10b8SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 268421e10b8SBarry Smith lastcol = col; 269421e10b8SBarry Smith while (high-low > 7) { 270421e10b8SBarry Smith t = (low+high)/2; 271421e10b8SBarry Smith if (rp[t] > bcol) high = t; 272421e10b8SBarry Smith else low = t; 273421e10b8SBarry Smith } 274421e10b8SBarry Smith for (i=low; i<high; i++) { 275421e10b8SBarry Smith if (rp[i] > bcol) break; 276421e10b8SBarry Smith if (rp[i] == bcol) { 277ccb205f8SBarry Smith /* printf("row %d col %d found i %d\n",brow,bcol,i);*/ 278421e10b8SBarry Smith goto noinsert1; 279421e10b8SBarry Smith } 280421e10b8SBarry Smith } 281421e10b8SBarry Smith if (nonew == 1) goto noinsert1; 282421e10b8SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 283421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat); 284421e10b8SBarry Smith N = nrow++ - 1; high++; 285421e10b8SBarry Smith /* shift up all the later entries in this row */ 286421e10b8SBarry Smith for (ii=N; ii>=i; ii--) { 287*a0e1a203SBarry Smith /* printf("shiffting N %d i %d ii %d brow %d bcol %d\n",N,i,ii,brow,bcol);*/ 288421e10b8SBarry Smith rp[ii+1] = rp[ii]; 289421e10b8SBarry Smith ap[ii+1] = ap[ii]; 290421e10b8SBarry Smith } 291421e10b8SBarry Smith if (N>=i) ap[i] = 0; 292421e10b8SBarry Smith rp[i] = bcol; 293421e10b8SBarry Smith a->nz++; 294421e10b8SBarry Smith noinsert1:; 295421e10b8SBarry Smith if (!*(ap+i)) { 296ccb205f8SBarry Smith /*printf("create matrix at i %d rw %d col %d\n",i,brow,bcol);*/ 297b0223207SBarry Smith if (A->symmetric && brow == bcol) { 2986e63c7a1SBarry Smith /* don't use SBAIJ since want to reorder in sparse factorization */ 2996e63c7a1SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr); 300b0223207SBarry Smith } else { 301421e10b8SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr); 302421e10b8SBarry Smith } 303b0223207SBarry Smith } 304ccb205f8SBarry Smith /* printf("numerical value at i %d row %d col %d cidx %d ridx %d value %g\n",i,brow,bcol,cidx,ridx,value);*/ 305421e10b8SBarry Smith ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr); 306421e10b8SBarry Smith low = i; 307421e10b8SBarry Smith } 308ccb205f8SBarry Smith /* printf("nrow for row %d %d\n",nrow,brow);*/ 309421e10b8SBarry Smith ailen[brow] = nrow; 310421e10b8SBarry Smith } 311*a0e1a203SBarry Smith /*printf("nz %d\n",a->nz);*/ 312421e10b8SBarry Smith A->same_nonzero = PETSC_FALSE; 313421e10b8SBarry Smith PetscFunctionReturn(0); 314421e10b8SBarry Smith } 315421e10b8SBarry Smith 316421e10b8SBarry Smith #undef __FUNCT__ 317421e10b8SBarry Smith #define __FUNCT__ "MatLoad_BlockMat" 318421e10b8SBarry Smith PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, MatType type,Mat *A) 319421e10b8SBarry Smith { 320421e10b8SBarry Smith PetscErrorCode ierr; 321421e10b8SBarry Smith Mat tmpA; 322*a0e1a203SBarry Smith PetscInt i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0; 323421e10b8SBarry Smith const PetscInt *cols; 324421e10b8SBarry Smith const PetscScalar *values; 3258cdf2d9bSBarry Smith PetscTruth flg,notdone; 3268cdf2d9bSBarry Smith Mat_SeqAIJ *a; 327*a0e1a203SBarry Smith Mat_BlockMat *amat; 328421e10b8SBarry Smith 329421e10b8SBarry Smith PetscFunctionBegin; 330421e10b8SBarry Smith ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr); 331421e10b8SBarry Smith 332421e10b8SBarry Smith ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr); 33377925062SSatish Balay ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr); 334421e10b8SBarry Smith ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr); 335290bbb0aSBarry Smith ierr = PetscOptionsName("-matload_symmetric","Store the matrix as symmetric","MatLoad",&flg);CHKERRQ(ierr); 3368cdf2d9bSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 3378cdf2d9bSBarry Smith 338*a0e1a203SBarry Smith /* Determine number of nonzero blocks for each block row */ 3398cdf2d9bSBarry Smith a = (Mat_SeqAIJ*) tmpA->data; 3408cdf2d9bSBarry Smith mbs = m/bs; 3418cdf2d9bSBarry Smith ierr = PetscMalloc3(mbs,PetscInt,&lens,bs,PetscInt,&ii,bs,PetscInt,&ilens);CHKERRQ(ierr); 3428cdf2d9bSBarry Smith ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr); 3438cdf2d9bSBarry Smith 3448cdf2d9bSBarry Smith for (i=0; i<mbs; i++) { 3458cdf2d9bSBarry Smith for (j=0; j<bs; j++) { 3468cdf2d9bSBarry Smith ii[j] = a->j + a->i[i*bs + j]; 3478cdf2d9bSBarry Smith ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 3488cdf2d9bSBarry Smith /* printf("j %d length %d\n",j,ilens[j]);*/ 3498cdf2d9bSBarry Smith } 350*a0e1a203SBarry Smith 351*a0e1a203SBarry Smith currentcol = -1; 3528cdf2d9bSBarry Smith notdone = PETSC_TRUE; 3538cdf2d9bSBarry Smith while (PETSC_TRUE) { 3548cdf2d9bSBarry Smith notdone = PETSC_FALSE; 3558cdf2d9bSBarry Smith nextcol = 1000000000; 3568cdf2d9bSBarry Smith for (j=0; j<bs; j++) { 3578cdf2d9bSBarry Smith /* printf("loop j %d length %d\n",j,ilens[j]); */ 358*a0e1a203SBarry Smith while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { 3598cdf2d9bSBarry Smith ii[j]++; 3608cdf2d9bSBarry Smith ilens[j]--; 3618cdf2d9bSBarry Smith } 3628cdf2d9bSBarry Smith if (ilens[j] > 0) { 3638cdf2d9bSBarry Smith notdone = PETSC_TRUE; 3648cdf2d9bSBarry Smith nextcol = PetscMin(nextcol,ii[j][0]/bs); 3658cdf2d9bSBarry Smith } 3668cdf2d9bSBarry Smith } 3678cdf2d9bSBarry Smith if (!notdone) break; 368*a0e1a203SBarry Smith printf("i %d currentcol %d lens[i] %d\n",i,currentcol,lens[i]); 369*a0e1a203SBarry Smith if (!flg || (nextcol >= i)) lens[i]++; 3708cdf2d9bSBarry Smith currentcol = nextcol; 3718cdf2d9bSBarry Smith } 372*a0e1a203SBarry Smith printf("len of i %d %d\n",i,lens[i]); 3738cdf2d9bSBarry Smith } 3748cdf2d9bSBarry Smith 3758cdf2d9bSBarry Smith ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,0,lens,A);CHKERRQ(ierr); 376290bbb0aSBarry Smith if (flg) { 377290bbb0aSBarry Smith ierr = MatSetOption(*A,MAT_SYMMETRIC);CHKERRQ(ierr); 378290bbb0aSBarry Smith } 379*a0e1a203SBarry Smith amat = (Mat_BlockMat*)(*A)->data; 380*a0e1a203SBarry Smith 381*a0e1a203SBarry Smith /* preallocate the submatrices */ 382*a0e1a203SBarry Smith ierr = PetscMalloc(bs*sizeof(PetscInt),&llens);CHKERRQ(ierr); 383*a0e1a203SBarry Smith for (i=0; i<mbs; i++) { /* loops for block rows */ 384*a0e1a203SBarry Smith for (j=0; j<bs; j++) { 385*a0e1a203SBarry Smith ii[j] = a->j + a->i[i*bs + j]; 386*a0e1a203SBarry Smith ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 387*a0e1a203SBarry Smith /* printf("j %d length %d\n",j,ilens[j]);*/ 388*a0e1a203SBarry Smith } 389*a0e1a203SBarry Smith 390*a0e1a203SBarry Smith currentcol = 1000000000; 391*a0e1a203SBarry Smith for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */ 392*a0e1a203SBarry Smith if (ilens[j] > 0) { 393*a0e1a203SBarry Smith currentcol = PetscMin(currentcol,ii[j][0]/bs); 394*a0e1a203SBarry Smith } 395*a0e1a203SBarry Smith } 396*a0e1a203SBarry Smith 397*a0e1a203SBarry Smith notdone = PETSC_TRUE; 398*a0e1a203SBarry Smith while (PETSC_TRUE) { /* loops over blocks in block row */ 399*a0e1a203SBarry Smith 400*a0e1a203SBarry Smith notdone = PETSC_FALSE; 401*a0e1a203SBarry Smith nextcol = 1000000000; 402*a0e1a203SBarry Smith ierr = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr); 403*a0e1a203SBarry Smith for (j=0; j<bs; j++) { /* loop over rows in block */ 404*a0e1a203SBarry Smith /*printf("loop j %d length %d\n",j,ilens[j]); */ 405*a0e1a203SBarry Smith /*printf("current col %d %d\n",currentcol,ii[j][0]/bs);*/ 406*a0e1a203SBarry Smith while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */ 407*a0e1a203SBarry Smith /*printf("j %d in llens[j] %d\n",j,llens[j]);*/ 408*a0e1a203SBarry Smith ii[j]++; 409*a0e1a203SBarry Smith ilens[j]--; 410*a0e1a203SBarry Smith llens[j]++; 411*a0e1a203SBarry Smith } 412*a0e1a203SBarry Smith if (ilens[j] > 0) { 413*a0e1a203SBarry Smith notdone = PETSC_TRUE; 414*a0e1a203SBarry Smith nextcol = PetscMin(nextcol,ii[j][0]/bs); 415*a0e1a203SBarry Smith } 416*a0e1a203SBarry Smith } 417*a0e1a203SBarry Smith printf("cnt %d llens %d %d %d %d\n",cnt,llens[0],llens[1],llens[2],llens[3]); 418*a0e1a203SBarry Smith if (cnt >= amat->maxnz) SETERRQ1(PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt); 419*a0e1a203SBarry Smith if (!flg || currentcol >= i) { 420*a0e1a203SBarry Smith amat->j[cnt] = currentcol; 421*a0e1a203SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr); 422*a0e1a203SBarry Smith /*printf("mat %p row %d col %d\n",amat->a[cnt-1],i,currentcol);*/ 423*a0e1a203SBarry Smith } 424*a0e1a203SBarry Smith 425*a0e1a203SBarry Smith if (!notdone) break; 426*a0e1a203SBarry Smith currentcol = nextcol; 427*a0e1a203SBarry Smith } 428*a0e1a203SBarry Smith amat->ilen[i] = lens[i]; 429*a0e1a203SBarry Smith } 430*a0e1a203SBarry Smith CHKMEMQ; 431*a0e1a203SBarry Smith /*printf("total cnt %d\n",cnt);*/ 432*a0e1a203SBarry Smith 433*a0e1a203SBarry Smith ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr); 434*a0e1a203SBarry Smith ierr = PetscFree(llens);CHKERRQ(ierr); 435*a0e1a203SBarry Smith 436*a0e1a203SBarry Smith /* copy over the matrix, one row at a time */ 437421e10b8SBarry Smith for (i=0; i<m; i++) { 438421e10b8SBarry Smith ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 439421e10b8SBarry Smith ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr); 440ccb205f8SBarry Smith ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 441421e10b8SBarry Smith } 442421e10b8SBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 443421e10b8SBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 444421e10b8SBarry Smith PetscFunctionReturn(0); 445421e10b8SBarry Smith } 446421e10b8SBarry Smith 447ccb205f8SBarry Smith #undef __FUNCT__ 448ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat" 449ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer) 450ccb205f8SBarry Smith { 451ccb205f8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 452ccb205f8SBarry Smith PetscErrorCode ierr; 453ccb205f8SBarry Smith const char *name; 454ccb205f8SBarry Smith PetscViewerFormat format; 455ccb205f8SBarry Smith 456ccb205f8SBarry Smith PetscFunctionBegin; 457ccb205f8SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 458ccb205f8SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 459ccb205f8SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 460ccb205f8SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr); 461ccb205f8SBarry Smith } 462ccb205f8SBarry Smith PetscFunctionReturn(0); 463ccb205f8SBarry Smith } 464421e10b8SBarry Smith 465421e10b8SBarry Smith #undef __FUNCT__ 466421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat" 467421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat) 468421e10b8SBarry Smith { 469421e10b8SBarry Smith PetscErrorCode ierr; 470421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 471ccb205f8SBarry Smith PetscInt i; 472421e10b8SBarry Smith 473421e10b8SBarry Smith PetscFunctionBegin; 474421e10b8SBarry Smith if (bmat->right) { 475421e10b8SBarry Smith ierr = VecDestroy(bmat->right);CHKERRQ(ierr); 476421e10b8SBarry Smith } 477421e10b8SBarry Smith if (bmat->left) { 478421e10b8SBarry Smith ierr = VecDestroy(bmat->left);CHKERRQ(ierr); 479421e10b8SBarry Smith } 480290bbb0aSBarry Smith if (bmat->middle) { 481290bbb0aSBarry Smith ierr = VecDestroy(bmat->middle);CHKERRQ(ierr); 482290bbb0aSBarry Smith } 4836e63c7a1SBarry Smith if (bmat->workb) { 4846e63c7a1SBarry Smith ierr = VecDestroy(bmat->workb);CHKERRQ(ierr); 4856e63c7a1SBarry Smith } 486ccb205f8SBarry Smith if (bmat->diags) { 487ccb205f8SBarry Smith for (i=0; i<mat->rmap.n/mat->rmap.bs; i++) { 488ccb205f8SBarry Smith if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);} 489ccb205f8SBarry Smith } 490ccb205f8SBarry Smith } 491ccb205f8SBarry Smith if (bmat->a) { 492ccb205f8SBarry Smith for (i=0; i<bmat->nz; i++) { 493ccb205f8SBarry Smith if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);} 494ccb205f8SBarry Smith } 495ccb205f8SBarry Smith } 496ccb205f8SBarry Smith ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 497421e10b8SBarry Smith ierr = PetscFree(bmat);CHKERRQ(ierr); 498421e10b8SBarry Smith PetscFunctionReturn(0); 499421e10b8SBarry Smith } 500421e10b8SBarry Smith 501421e10b8SBarry Smith #undef __FUNCT__ 502421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat" 503421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 504421e10b8SBarry Smith { 505421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 506421e10b8SBarry Smith PetscErrorCode ierr; 507421e10b8SBarry Smith PetscScalar *xx,*yy; 508421e10b8SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j; 509421e10b8SBarry Smith Mat *aa; 510421e10b8SBarry Smith 511421e10b8SBarry Smith PetscFunctionBegin; 512ccb205f8SBarry Smith CHKMEMQ; 513421e10b8SBarry Smith /* 514421e10b8SBarry Smith Standard CSR multiply except each entry is a Mat 515421e10b8SBarry Smith */ 516421e10b8SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 517ccb205f8SBarry Smith 518421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 519421e10b8SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 520421e10b8SBarry Smith aj = bmat->j; 521421e10b8SBarry Smith aa = bmat->a; 522421e10b8SBarry Smith ii = bmat->i; 523421e10b8SBarry Smith for (i=0; i<m; i++) { 524421e10b8SBarry Smith jrow = ii[i]; 525ccb205f8SBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 526421e10b8SBarry Smith n = ii[i+1] - jrow; 527421e10b8SBarry Smith for (j=0; j<n; j++) { 528421e10b8SBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 529ccb205f8SBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 530421e10b8SBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 531421e10b8SBarry Smith jrow++; 532421e10b8SBarry Smith } 533421e10b8SBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 534421e10b8SBarry Smith } 535421e10b8SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 536421e10b8SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 537ccb205f8SBarry Smith CHKMEMQ; 538421e10b8SBarry Smith PetscFunctionReturn(0); 539421e10b8SBarry Smith } 540421e10b8SBarry Smith 541421e10b8SBarry Smith #undef __FUNCT__ 542290bbb0aSBarry Smith #define __FUNCT__ "MatMult_BlockMat_Symmetric" 543290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y) 544290bbb0aSBarry Smith { 545290bbb0aSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 546290bbb0aSBarry Smith PetscErrorCode ierr; 547290bbb0aSBarry Smith PetscScalar *xx,*yy; 548290bbb0aSBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j; 549290bbb0aSBarry Smith Mat *aa; 550290bbb0aSBarry Smith 551290bbb0aSBarry Smith PetscFunctionBegin; 552290bbb0aSBarry Smith CHKMEMQ; 553290bbb0aSBarry Smith /* 554290bbb0aSBarry Smith Standard CSR multiply except each entry is a Mat 555290bbb0aSBarry Smith */ 556290bbb0aSBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 557290bbb0aSBarry Smith 558290bbb0aSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 559290bbb0aSBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 560290bbb0aSBarry Smith aj = bmat->j; 561290bbb0aSBarry Smith aa = bmat->a; 562290bbb0aSBarry Smith ii = bmat->i; 563290bbb0aSBarry Smith for (i=0; i<m; i++) { 564290bbb0aSBarry Smith jrow = ii[i]; 565290bbb0aSBarry Smith n = ii[i+1] - jrow; 566290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 567290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr); 568290bbb0aSBarry Smith /* if we ALWAYS required a diagonal entry then could remove this if test */ 569290bbb0aSBarry Smith if (aj[jrow] == i) { 570290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 571290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 572290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 573290bbb0aSBarry Smith jrow++; 574290bbb0aSBarry Smith n--; 575290bbb0aSBarry Smith } 576290bbb0aSBarry Smith for (j=0; j<n; j++) { 577290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); /* upper triangular part */ 578290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 579290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 580290bbb0aSBarry Smith 581290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr); /* lower triangular part */ 582290bbb0aSBarry Smith ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr); 583290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 584290bbb0aSBarry Smith jrow++; 585290bbb0aSBarry Smith } 586290bbb0aSBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 587290bbb0aSBarry Smith ierr = VecResetArray(bmat->middle);CHKERRQ(ierr); 588290bbb0aSBarry Smith } 589290bbb0aSBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 590290bbb0aSBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 591290bbb0aSBarry Smith CHKMEMQ; 592290bbb0aSBarry Smith PetscFunctionReturn(0); 593290bbb0aSBarry Smith } 594290bbb0aSBarry Smith 595290bbb0aSBarry Smith #undef __FUNCT__ 596421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat" 597421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 598421e10b8SBarry Smith { 599421e10b8SBarry Smith PetscFunctionBegin; 600421e10b8SBarry Smith PetscFunctionReturn(0); 601421e10b8SBarry Smith } 602421e10b8SBarry Smith 603421e10b8SBarry Smith #undef __FUNCT__ 604421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat" 605421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 606421e10b8SBarry Smith { 607421e10b8SBarry Smith PetscFunctionBegin; 608421e10b8SBarry Smith PetscFunctionReturn(0); 609421e10b8SBarry Smith } 610421e10b8SBarry Smith 611421e10b8SBarry Smith #undef __FUNCT__ 612421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat" 613421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 614421e10b8SBarry Smith { 615421e10b8SBarry Smith PetscFunctionBegin; 616421e10b8SBarry Smith PetscFunctionReturn(0); 617421e10b8SBarry Smith } 618421e10b8SBarry Smith 619ccb205f8SBarry Smith /* 620ccb205f8SBarry Smith Adds diagonal pointers to sparse matrix structure. 621ccb205f8SBarry Smith */ 622ccb205f8SBarry Smith #undef __FUNCT__ 623ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat" 624ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) 625ccb205f8SBarry Smith { 626ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 627ccb205f8SBarry Smith PetscErrorCode ierr; 628ccb205f8SBarry Smith PetscInt i,j,mbs = A->rmap.n/A->rmap.bs; 629ccb205f8SBarry Smith 630ccb205f8SBarry Smith PetscFunctionBegin; 631ccb205f8SBarry Smith if (!a->diag) { 632ccb205f8SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 633ccb205f8SBarry Smith } 634ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 635ccb205f8SBarry Smith a->diag[i] = a->i[i+1]; 636ccb205f8SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 637ccb205f8SBarry Smith if (a->j[j] == i) { 638ccb205f8SBarry Smith a->diag[i] = j; 639ccb205f8SBarry Smith break; 640ccb205f8SBarry Smith } 641ccb205f8SBarry Smith } 642ccb205f8SBarry Smith } 643ccb205f8SBarry Smith PetscFunctionReturn(0); 644ccb205f8SBarry Smith } 645ccb205f8SBarry Smith 646ccb205f8SBarry Smith #undef __FUNCT__ 647ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat" 648ccb205f8SBarry Smith PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B) 649ccb205f8SBarry Smith { 650ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 651ccb205f8SBarry Smith Mat_SeqAIJ *c; 652ccb205f8SBarry Smith PetscErrorCode ierr; 653ccb205f8SBarry Smith PetscInt i,k,first,step,lensi,nrows,ncols; 654ccb205f8SBarry Smith PetscInt *j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen; 655ccb205f8SBarry Smith PetscScalar *a_new,value; 656ccb205f8SBarry Smith Mat C,*aa = a->a; 657ccb205f8SBarry Smith PetscTruth stride,equal; 658ccb205f8SBarry Smith 659ccb205f8SBarry Smith PetscFunctionBegin; 660ccb205f8SBarry Smith ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr); 661ccb205f8SBarry Smith if (!equal) SETERRQ(PETSC_ERR_SUP,"Only for idential column and row indices"); 662ccb205f8SBarry Smith ierr = ISStride(iscol,&stride);CHKERRQ(ierr); 663ccb205f8SBarry Smith if (!stride) SETERRQ(PETSC_ERR_SUP,"Only for stride indices"); 664ccb205f8SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 665ccb205f8SBarry Smith if (step != A->rmap.bs) SETERRQ(PETSC_ERR_SUP,"Can only select one entry from each block"); 666ccb205f8SBarry Smith 667ccb205f8SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 668ccb205f8SBarry Smith ncols = nrows; 669ccb205f8SBarry Smith 670ccb205f8SBarry Smith /* create submatrix */ 671ccb205f8SBarry Smith if (scall == MAT_REUSE_MATRIX) { 672ccb205f8SBarry Smith PetscInt n_cols,n_rows; 673ccb205f8SBarry Smith C = *B; 674ccb205f8SBarry Smith ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr); 675ccb205f8SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 676ccb205f8SBarry Smith ierr = MatZeroEntries(C);CHKERRQ(ierr); 677ccb205f8SBarry Smith } else { 678ccb205f8SBarry Smith ierr = MatCreate(A->comm,&C);CHKERRQ(ierr); 679ccb205f8SBarry Smith ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 6806e63c7a1SBarry Smith if (A->symmetric) { 6816e63c7a1SBarry Smith ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 6826e63c7a1SBarry Smith } else { 683ccb205f8SBarry Smith ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr); 6846e63c7a1SBarry Smith } 6856e63c7a1SBarry Smith ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr); 6866e63c7a1SBarry Smith ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr); 687ccb205f8SBarry Smith } 688ccb205f8SBarry Smith c = (Mat_SeqAIJ*)C->data; 689ccb205f8SBarry Smith 690ccb205f8SBarry Smith /* loop over rows inserting into submatrix */ 691ccb205f8SBarry Smith a_new = c->a; 692ccb205f8SBarry Smith j_new = c->j; 693ccb205f8SBarry Smith i_new = c->i; 694ccb205f8SBarry Smith 695ccb205f8SBarry Smith for (i=0; i<nrows; i++) { 696ccb205f8SBarry Smith ii = ai[i]; 697ccb205f8SBarry Smith lensi = ailen[i]; 698ccb205f8SBarry Smith for (k=0; k<lensi; k++) { 699ccb205f8SBarry Smith *j_new++ = *aj++; 700ccb205f8SBarry Smith ierr = MatGetValue(*aa++,first,first,value);CHKERRQ(ierr); 701ccb205f8SBarry Smith *a_new++ = value; 702ccb205f8SBarry Smith } 703ccb205f8SBarry Smith i_new[i+1] = i_new[i] + lensi; 704ccb205f8SBarry Smith c->ilen[i] = lensi; 705ccb205f8SBarry Smith } 706ccb205f8SBarry Smith 707ccb205f8SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 708ccb205f8SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 709ccb205f8SBarry Smith *B = C; 710ccb205f8SBarry Smith PetscFunctionReturn(0); 711ccb205f8SBarry Smith } 712ccb205f8SBarry Smith 713421e10b8SBarry Smith #undef __FUNCT__ 714421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat" 715421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 716421e10b8SBarry Smith { 717421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 718421e10b8SBarry Smith PetscErrorCode ierr; 719421e10b8SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 720ccb205f8SBarry Smith PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0; 721421e10b8SBarry Smith Mat *aa = a->a,*ap; 722421e10b8SBarry Smith 723421e10b8SBarry Smith PetscFunctionBegin; 724421e10b8SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 725421e10b8SBarry Smith 726421e10b8SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 727421e10b8SBarry Smith for (i=1; i<m; i++) { 728421e10b8SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 729421e10b8SBarry Smith fshift += imax[i-1] - ailen[i-1]; 730421e10b8SBarry Smith rmax = PetscMax(rmax,ailen[i]); 731421e10b8SBarry Smith if (fshift) { 732421e10b8SBarry Smith ip = aj + ai[i] ; 733421e10b8SBarry Smith ap = aa + ai[i] ; 734421e10b8SBarry Smith N = ailen[i]; 735421e10b8SBarry Smith for (j=0; j<N; j++) { 736421e10b8SBarry Smith ip[j-fshift] = ip[j]; 737421e10b8SBarry Smith ap[j-fshift] = ap[j]; 738421e10b8SBarry Smith } 739421e10b8SBarry Smith } 740421e10b8SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 741421e10b8SBarry Smith } 742421e10b8SBarry Smith if (m) { 743421e10b8SBarry Smith fshift += imax[m-1] - ailen[m-1]; 744421e10b8SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 745421e10b8SBarry Smith } 746421e10b8SBarry Smith /* reset ilen and imax for each row */ 747421e10b8SBarry Smith for (i=0; i<m; i++) { 748421e10b8SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 749421e10b8SBarry Smith } 750421e10b8SBarry Smith a->nz = ai[m]; 751ccb205f8SBarry Smith for (i=0; i<a->nz; i++) { 752ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG) 753ccb205f8SBarry Smith if (!aa[i]) SETERRQ3(PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz); 754ccb205f8SBarry Smith #endif 755*a0e1a203SBarry Smith /*printf("mat assembly %p col %d\n",aa[i],aj[i]);*/ 756ccb205f8SBarry Smith ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 757ccb205f8SBarry Smith ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 758ccb205f8SBarry Smith } 759ccb205f8SBarry Smith CHKMEMQ; 7608cdf2d9bSBarry Smith ierr = PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n/A->cmap.bs,fshift,a->nz);CHKERRQ(ierr); 761421e10b8SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); 762421e10b8SBarry Smith ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); 763421e10b8SBarry Smith a->reallocs = 0; 764421e10b8SBarry Smith A->info.nz_unneeded = (double)fshift; 765421e10b8SBarry Smith a->rmax = rmax; 766421e10b8SBarry Smith 767421e10b8SBarry Smith A->same_nonzero = PETSC_TRUE; 768ccb205f8SBarry Smith ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr); 769421e10b8SBarry Smith PetscFunctionReturn(0); 770421e10b8SBarry Smith } 771421e10b8SBarry Smith 772290bbb0aSBarry Smith #undef __FUNCT__ 773290bbb0aSBarry Smith #define __FUNCT__ "MatSetOption_BlockMat" 774290bbb0aSBarry Smith PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt) 775290bbb0aSBarry Smith { 776290bbb0aSBarry Smith PetscFunctionBegin; 777290bbb0aSBarry Smith if (opt == MAT_SYMMETRIC) { 778290bbb0aSBarry Smith A->ops->relax = MatRelax_BlockMat_Symmetric; 779290bbb0aSBarry Smith A->ops->mult = MatMult_BlockMat_Symmetric; 780290bbb0aSBarry Smith } else { 781290bbb0aSBarry Smith PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]); 782290bbb0aSBarry Smith } 783290bbb0aSBarry Smith PetscFunctionReturn(0); 784290bbb0aSBarry Smith } 785290bbb0aSBarry Smith 786290bbb0aSBarry Smith 787421e10b8SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 788421e10b8SBarry Smith 0, 789421e10b8SBarry Smith 0, 790421e10b8SBarry Smith MatMult_BlockMat, 791421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat, 792421e10b8SBarry Smith MatMultTranspose_BlockMat, 793421e10b8SBarry Smith MatMultTransposeAdd_BlockMat, 794421e10b8SBarry Smith 0, 795421e10b8SBarry Smith 0, 796421e10b8SBarry Smith 0, 797421e10b8SBarry Smith /*10*/ 0, 798421e10b8SBarry Smith 0, 799421e10b8SBarry Smith 0, 800ccb205f8SBarry Smith MatRelax_BlockMat, 801421e10b8SBarry Smith 0, 802421e10b8SBarry Smith /*15*/ 0, 803421e10b8SBarry Smith 0, 804421e10b8SBarry Smith 0, 805421e10b8SBarry Smith 0, 806421e10b8SBarry Smith 0, 807421e10b8SBarry Smith /*20*/ 0, 808421e10b8SBarry Smith MatAssemblyEnd_BlockMat, 809421e10b8SBarry Smith 0, 810290bbb0aSBarry Smith MatSetOption_BlockMat, 811421e10b8SBarry Smith 0, 812421e10b8SBarry Smith /*25*/ 0, 813421e10b8SBarry Smith 0, 814421e10b8SBarry Smith 0, 815421e10b8SBarry Smith 0, 816421e10b8SBarry Smith 0, 817421e10b8SBarry Smith /*30*/ 0, 818421e10b8SBarry Smith 0, 819421e10b8SBarry Smith 0, 820421e10b8SBarry Smith 0, 821421e10b8SBarry Smith 0, 822421e10b8SBarry Smith /*35*/ 0, 823421e10b8SBarry Smith 0, 824421e10b8SBarry Smith 0, 825421e10b8SBarry Smith 0, 826421e10b8SBarry Smith 0, 827421e10b8SBarry Smith /*40*/ 0, 828421e10b8SBarry Smith 0, 829421e10b8SBarry Smith 0, 830421e10b8SBarry Smith 0, 831421e10b8SBarry Smith 0, 832421e10b8SBarry Smith /*45*/ 0, 833421e10b8SBarry Smith 0, 834421e10b8SBarry Smith 0, 835421e10b8SBarry Smith 0, 836421e10b8SBarry Smith 0, 8378cdf2d9bSBarry Smith /*50*/ 0, 838421e10b8SBarry Smith 0, 839421e10b8SBarry Smith 0, 840421e10b8SBarry Smith 0, 841421e10b8SBarry Smith 0, 842421e10b8SBarry Smith /*55*/ 0, 843421e10b8SBarry Smith 0, 844421e10b8SBarry Smith 0, 845421e10b8SBarry Smith 0, 846421e10b8SBarry Smith 0, 847ccb205f8SBarry Smith /*60*/ MatGetSubMatrix_BlockMat, 848421e10b8SBarry Smith MatDestroy_BlockMat, 849ccb205f8SBarry Smith MatView_BlockMat, 850421e10b8SBarry Smith 0, 851421e10b8SBarry Smith 0, 852421e10b8SBarry Smith /*65*/ 0, 853421e10b8SBarry Smith 0, 854421e10b8SBarry Smith 0, 855421e10b8SBarry Smith 0, 856421e10b8SBarry Smith 0, 857421e10b8SBarry Smith /*70*/ 0, 858421e10b8SBarry Smith 0, 859421e10b8SBarry Smith 0, 860421e10b8SBarry Smith 0, 861421e10b8SBarry Smith 0, 862421e10b8SBarry Smith /*75*/ 0, 863421e10b8SBarry Smith 0, 864421e10b8SBarry Smith 0, 865421e10b8SBarry Smith 0, 866421e10b8SBarry Smith 0, 867421e10b8SBarry Smith /*80*/ 0, 868421e10b8SBarry Smith 0, 869421e10b8SBarry Smith 0, 870421e10b8SBarry Smith 0, 871421e10b8SBarry Smith MatLoad_BlockMat, 872421e10b8SBarry Smith /*85*/ 0, 873421e10b8SBarry Smith 0, 874421e10b8SBarry Smith 0, 875421e10b8SBarry Smith 0, 876421e10b8SBarry Smith 0, 877421e10b8SBarry Smith /*90*/ 0, 878421e10b8SBarry Smith 0, 879421e10b8SBarry Smith 0, 880421e10b8SBarry Smith 0, 881421e10b8SBarry Smith 0, 882421e10b8SBarry Smith /*95*/ 0, 883421e10b8SBarry Smith 0, 884421e10b8SBarry Smith 0, 885421e10b8SBarry Smith 0}; 886421e10b8SBarry Smith 8878cdf2d9bSBarry Smith #undef __FUNCT__ 8888cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation" 8898cdf2d9bSBarry Smith /*@C 8908cdf2d9bSBarry Smith MatBlockMatSetPreallocation - For good matrix assembly performance 8918cdf2d9bSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 8928cdf2d9bSBarry Smith (or the array nnz). By setting these parameters accurately, performance 8938cdf2d9bSBarry Smith during matrix assembly can be increased by more than a factor of 50. 8948cdf2d9bSBarry Smith 8958cdf2d9bSBarry Smith Collective on MPI_Comm 8968cdf2d9bSBarry Smith 8978cdf2d9bSBarry Smith Input Parameters: 8988cdf2d9bSBarry Smith + B - The matrix 8998cdf2d9bSBarry Smith . bs - size of each block in matrix 9008cdf2d9bSBarry Smith . nz - number of nonzeros per block row (same for all rows) 9018cdf2d9bSBarry Smith - nnz - array containing the number of nonzeros in the various block rows 9028cdf2d9bSBarry Smith (possibly different for each row) or PETSC_NULL 9038cdf2d9bSBarry Smith 9048cdf2d9bSBarry Smith Notes: 9058cdf2d9bSBarry Smith If nnz is given then nz is ignored 9068cdf2d9bSBarry Smith 9078cdf2d9bSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 9088cdf2d9bSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 9098cdf2d9bSBarry Smith allocation. For large problems you MUST preallocate memory or you 9108cdf2d9bSBarry Smith will get TERRIBLE performance, see the users' manual chapter on matrices. 9118cdf2d9bSBarry Smith 9128cdf2d9bSBarry Smith Level: intermediate 9138cdf2d9bSBarry Smith 9148cdf2d9bSBarry Smith .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues() 9158cdf2d9bSBarry Smith 9168cdf2d9bSBarry Smith @*/ 9178cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 9188cdf2d9bSBarry Smith { 9198cdf2d9bSBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 9208cdf2d9bSBarry Smith 9218cdf2d9bSBarry Smith PetscFunctionBegin; 9228cdf2d9bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)B,"MatBlockMatSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 9238cdf2d9bSBarry Smith if (f) { 9248cdf2d9bSBarry Smith ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 9258cdf2d9bSBarry Smith } 9268cdf2d9bSBarry Smith PetscFunctionReturn(0); 9278cdf2d9bSBarry Smith } 9288cdf2d9bSBarry Smith 9298cdf2d9bSBarry Smith EXTERN_C_BEGIN 9308cdf2d9bSBarry Smith #undef __FUNCT__ 9318cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation_BlockMat" 9328cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz) 9338cdf2d9bSBarry Smith { 9348cdf2d9bSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 9358cdf2d9bSBarry Smith PetscErrorCode ierr; 9368cdf2d9bSBarry Smith PetscInt i; 9378cdf2d9bSBarry Smith 9388cdf2d9bSBarry Smith PetscFunctionBegin; 9398cdf2d9bSBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Block size given %D must be great than zero",bs); 9408cdf2d9bSBarry Smith if (A->rmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap.n); 9418cdf2d9bSBarry Smith if (A->cmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap.n); 9428cdf2d9bSBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 9438cdf2d9bSBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 9448cdf2d9bSBarry Smith if (nnz) { 9458cdf2d9bSBarry Smith for (i=0; i<A->rmap.n/bs; i++) { 9468cdf2d9bSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 9478cdf2d9bSBarry Smith if (nnz[i] > A->cmap.n/bs) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],A->cmap.n/bs); 9488cdf2d9bSBarry Smith } 9498cdf2d9bSBarry Smith } 9508cdf2d9bSBarry Smith A->rmap.bs = A->cmap.bs = bs; 9518cdf2d9bSBarry Smith bmat->mbs = A->rmap.n/bs; 9528cdf2d9bSBarry Smith 953*a0e1a203SBarry Smith /*printf("A->rmap.n%d %d\n",A->rmap.n, bs);*/ 9548cdf2d9bSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr); 9558cdf2d9bSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->middle);CHKERRQ(ierr); 9568cdf2d9bSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr); 9578cdf2d9bSBarry Smith 9588cdf2d9bSBarry Smith ierr = PetscMalloc2(A->rmap.n,PetscInt,&bmat->imax,A->rmap.n,PetscInt,&bmat->ilen);CHKERRQ(ierr); 9598cdf2d9bSBarry Smith if (nnz) { 9608cdf2d9bSBarry Smith nz = 0; 9618cdf2d9bSBarry Smith for (i=0; i<A->rmap.n/A->rmap.bs; i++) { 9628cdf2d9bSBarry Smith bmat->imax[i] = nnz[i]; 9638cdf2d9bSBarry Smith nz += nnz[i]; 9648cdf2d9bSBarry Smith } 9658cdf2d9bSBarry Smith } else { 9668cdf2d9bSBarry Smith SETERRQ(PETSC_ERR_SUP,"Currently requires block row by row preallocation"); 9678cdf2d9bSBarry Smith } 9688cdf2d9bSBarry Smith 969*a0e1a203SBarry Smith /*printf("nz in p;real %d\n",nz);*/ 970*a0e1a203SBarry Smith 9718cdf2d9bSBarry Smith /* bmat->ilen will count nonzeros in each row so far. */ 9728cdf2d9bSBarry Smith for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;} 9738cdf2d9bSBarry Smith 9748cdf2d9bSBarry Smith /* allocate the matrix space */ 9758cdf2d9bSBarry Smith ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap.n+1,PetscInt,&bmat->i);CHKERRQ(ierr); 9768cdf2d9bSBarry Smith bmat->i[0] = 0; 9778cdf2d9bSBarry Smith for (i=1; i<bmat->mbs+1; i++) { 9788cdf2d9bSBarry Smith bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1]; 9798cdf2d9bSBarry Smith } 9808cdf2d9bSBarry Smith bmat->singlemalloc = PETSC_TRUE; 9818cdf2d9bSBarry Smith bmat->free_a = PETSC_TRUE; 9828cdf2d9bSBarry Smith bmat->free_ij = PETSC_TRUE; 9838cdf2d9bSBarry Smith 9848cdf2d9bSBarry Smith bmat->nz = 0; 9858cdf2d9bSBarry Smith bmat->maxnz = nz; 9868cdf2d9bSBarry Smith A->info.nz_unneeded = (double)bmat->maxnz; 9878cdf2d9bSBarry Smith 9888cdf2d9bSBarry Smith PetscFunctionReturn(0); 9898cdf2d9bSBarry Smith } 9908cdf2d9bSBarry Smith EXTERN_C_END 9918cdf2d9bSBarry Smith 992421e10b8SBarry Smith /*MC 993421e10b8SBarry Smith MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 994421e10b8SBarry Smith consisting of (usually) sparse blocks. 995421e10b8SBarry Smith 996421e10b8SBarry Smith Level: advanced 997421e10b8SBarry Smith 998421e10b8SBarry Smith .seealso: MatCreateBlockMat() 999421e10b8SBarry Smith 1000421e10b8SBarry Smith M*/ 1001421e10b8SBarry Smith 1002421e10b8SBarry Smith EXTERN_C_BEGIN 1003421e10b8SBarry Smith #undef __FUNCT__ 1004421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat" 1005421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A) 1006421e10b8SBarry Smith { 1007421e10b8SBarry Smith Mat_BlockMat *b; 1008421e10b8SBarry Smith PetscErrorCode ierr; 1009421e10b8SBarry Smith 1010421e10b8SBarry Smith PetscFunctionBegin; 1011421e10b8SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1012421e10b8SBarry Smith ierr = PetscNew(Mat_BlockMat,&b);CHKERRQ(ierr); 1013421e10b8SBarry Smith 1014421e10b8SBarry Smith A->data = (void*)b; 1015421e10b8SBarry Smith 1016421e10b8SBarry Smith ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr); 1017421e10b8SBarry Smith ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr); 1018421e10b8SBarry Smith 1019421e10b8SBarry Smith A->assembled = PETSC_TRUE; 1020421e10b8SBarry Smith A->preallocated = PETSC_FALSE; 1021421e10b8SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr); 1022290bbb0aSBarry Smith 1023290bbb0aSBarry Smith ierr = PetscOptionsBegin(A->comm,A->prefix,"Matrix Option","Mat");CHKERRQ(ierr); 1024290bbb0aSBarry Smith ierr = PetscOptionsEnd(); 1025290bbb0aSBarry Smith 10268cdf2d9bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatBlockMatSetPreallocation_C", 10278cdf2d9bSBarry Smith "MatBlockMatSetPreallocation_BlockMat", 10288cdf2d9bSBarry Smith MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr); 10298cdf2d9bSBarry Smith 1030421e10b8SBarry Smith PetscFunctionReturn(0); 1031421e10b8SBarry Smith } 1032421e10b8SBarry Smith EXTERN_C_END 1033421e10b8SBarry Smith 1034421e10b8SBarry Smith #undef __FUNCT__ 1035421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat" 1036421e10b8SBarry Smith /*@C 1037421e10b8SBarry Smith MatCreateBlockMat - Creates a new matrix based sparse Mat storage 1038421e10b8SBarry Smith 1039421e10b8SBarry Smith Collective on MPI_Comm 1040421e10b8SBarry Smith 1041421e10b8SBarry Smith Input Parameters: 1042421e10b8SBarry Smith + comm - MPI communicator 1043421e10b8SBarry Smith . m - number of rows 1044421e10b8SBarry Smith . n - number of columns 10458cdf2d9bSBarry Smith . bs - size of each submatrix 10468cdf2d9bSBarry Smith . nz - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known) 10478cdf2d9bSBarry Smith - nnz - expected number of nonzers per block row if known (use PETSC_NULL otherwise) 1048421e10b8SBarry Smith 1049421e10b8SBarry Smith 1050421e10b8SBarry Smith Output Parameter: 1051421e10b8SBarry Smith . A - the matrix 1052421e10b8SBarry Smith 1053421e10b8SBarry Smith Level: intermediate 1054421e10b8SBarry Smith 1055421e10b8SBarry Smith PETSc requires that matrices and vectors being used for certain 1056421e10b8SBarry Smith operations are partitioned accordingly. For example, when 1057421e10b8SBarry Smith creating a bmat matrix, A, that supports parallel matrix-vector 1058421e10b8SBarry Smith products using MatMult(A,x,y) the user should set the number 1059421e10b8SBarry Smith of local matrix rows to be the number of local elements of the 1060421e10b8SBarry Smith corresponding result vector, y. Note that this is information is 1061421e10b8SBarry Smith required for use of the matrix interface routines, even though 1062421e10b8SBarry Smith the bmat matrix may not actually be physically partitioned. 1063421e10b8SBarry Smith For example, 1064421e10b8SBarry Smith 1065421e10b8SBarry Smith .keywords: matrix, bmat, create 1066421e10b8SBarry Smith 1067421e10b8SBarry Smith .seealso: MATBLOCKMAT 1068421e10b8SBarry Smith @*/ 10698cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A) 1070421e10b8SBarry Smith { 1071421e10b8SBarry Smith PetscErrorCode ierr; 1072421e10b8SBarry Smith 1073421e10b8SBarry Smith PetscFunctionBegin; 1074421e10b8SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1075421e10b8SBarry Smith ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1076421e10b8SBarry Smith ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr); 10778cdf2d9bSBarry Smith ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 1078421e10b8SBarry Smith PetscFunctionReturn(0); 1079421e10b8SBarry Smith } 1080421e10b8SBarry Smith 1081421e10b8SBarry Smith 1082421e10b8SBarry Smith 1083