1421e10b8SBarry Smith #define PETSCMAT_DLL 2421e10b8SBarry Smith 3421e10b8SBarry Smith /* 4421e10b8SBarry Smith This provides a matrix that consists of Mats 5421e10b8SBarry Smith */ 6421e10b8SBarry Smith 77c4f633dSBarry Smith #include "private/matimpl.h" /*I "petscmat.h" I*/ 87c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h" /* use the common AIJ data-structure */ 9ccb205f8SBarry Smith #include "petscksp.h" 10421e10b8SBarry Smith 11421e10b8SBarry Smith typedef struct { 12421e10b8SBarry Smith SEQAIJHEADER(Mat); 13421e10b8SBarry Smith SEQBAIJHEADER; 14ccb205f8SBarry Smith Mat *diags; 15421e10b8SBarry Smith 166e63c7a1SBarry Smith Vec left,right,middle,workb; /* dummy vectors to perform local parts of product */ 17421e10b8SBarry Smith } Mat_BlockMat; 18421e10b8SBarry Smith 19a5973382SShri Abhyankar EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt*); 20a5973382SShri Abhyankar 21421e10b8SBarry Smith #undef __FUNCT__ 2241f059aeSBarry Smith #define __FUNCT__ "MatSOR_BlockMat_Symmetric" 2341f059aeSBarry Smith PetscErrorCode MatSOR_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; 30d0f46423SBarry 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; 39e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 40e32f2f54SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 41e32f2f54SBarry Smith if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 42e32f2f54SBarry Smith if (fshift) SETERRQ(PETSC_COMM_SELF,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)) 44e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,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++) { 502692d6eeSBarry Smith ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr); 51719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info);CHKERRQ(ierr); 52719d5645SBarry Smith ierr = MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);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); 643649974fSBarry Smith ierr = VecGetArrayRead(a->workb,&b);CHKERRQ(ierr); 65290bbb0aSBarry Smith 6641f059aeSBarry Smith /* need to add code for when initial guess is zero, see MatSOR_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 8841f059aeSBarry Smith /* now adjust right hand side, see MatSOR_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); 1243649974fSBarry Smith ierr = VecRestoreArrayRead(a->workb,&b);CHKERRQ(ierr); 125290bbb0aSBarry Smith PetscFunctionReturn(0); 126290bbb0aSBarry Smith } 127290bbb0aSBarry Smith 128290bbb0aSBarry Smith #undef __FUNCT__ 12941f059aeSBarry Smith #define __FUNCT__ "MatSOR_BlockMat" 13041f059aeSBarry Smith PetscErrorCode MatSOR_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; 137d0f46423SBarry 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; 146e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 147e32f2f54SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 148e32f2f54SBarry Smith if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 149e32f2f54SBarry Smith if (fshift) SETERRQ(PETSC_COMM_SELF,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++) { 1552692d6eeSBarry Smith ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr); 156719d5645SBarry Smith ierr = MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info);CHKERRQ(ierr); 157719d5645SBarry Smith ierr = MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);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); 1663649974fSBarry Smith ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 167ccb205f8SBarry Smith 16841f059aeSBarry Smith /* need to add code for when initial guess is zero, see MatSOR_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); 2213649974fSBarry Smith ierr = VecRestoreArrayRead(bb,&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; 232d0f46423SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol; 233421e10b8SBarry Smith PetscErrorCode ierr; 234421e10b8SBarry Smith PetscInt ridx,cidx; 235ace3abfcSBarry Smith PetscBool roworiented=a->roworiented; 236421e10b8SBarry Smith MatScalar value; 237421e10b8SBarry Smith Mat *ap,*aa = a->a; 238421e10b8SBarry Smith 239421e10b8SBarry Smith PetscFunctionBegin; 24071fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 241421e10b8SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 242421e10b8SBarry Smith row = im[k]; 243421e10b8SBarry Smith brow = row/bs; 244421e10b8SBarry Smith if (row < 0) continue; 245421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 246e32f2f54SBarry Smith if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1); 247421e10b8SBarry Smith #endif 248421e10b8SBarry Smith rp = aj + ai[brow]; 249421e10b8SBarry Smith ap = aa + ai[brow]; 250421e10b8SBarry Smith rmax = imax[brow]; 251421e10b8SBarry Smith nrow = ailen[brow]; 252421e10b8SBarry Smith low = 0; 253421e10b8SBarry Smith high = nrow; 254421e10b8SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 255421e10b8SBarry Smith if (in[l] < 0) continue; 256421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 257e32f2f54SBarry Smith if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1); 258421e10b8SBarry Smith #endif 259421e10b8SBarry Smith col = in[l]; bcol = col/bs; 2606e63c7a1SBarry Smith if (A->symmetric && brow > bcol) continue; 261421e10b8SBarry Smith ridx = row % bs; cidx = col % bs; 262421e10b8SBarry Smith if (roworiented) { 263421e10b8SBarry Smith value = v[l + k*n]; 264421e10b8SBarry Smith } else { 265421e10b8SBarry Smith value = v[k + l*m]; 266421e10b8SBarry Smith } 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) { 277421e10b8SBarry Smith goto noinsert1; 278421e10b8SBarry Smith } 279421e10b8SBarry Smith } 280421e10b8SBarry Smith if (nonew == 1) goto noinsert1; 281e32f2f54SBarry Smith if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 282fef13f97SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat); 283421e10b8SBarry Smith N = nrow++ - 1; high++; 284421e10b8SBarry Smith /* shift up all the later entries in this row */ 285421e10b8SBarry Smith for (ii=N; ii>=i; ii--) { 286421e10b8SBarry Smith rp[ii+1] = rp[ii]; 287421e10b8SBarry Smith ap[ii+1] = ap[ii]; 288421e10b8SBarry Smith } 289421e10b8SBarry Smith if (N>=i) ap[i] = 0; 290421e10b8SBarry Smith rp[i] = bcol; 291421e10b8SBarry Smith a->nz++; 292421e10b8SBarry Smith noinsert1:; 293421e10b8SBarry Smith if (!*(ap+i)) { 2946e63c7a1SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr); 295b0223207SBarry Smith } 296421e10b8SBarry Smith ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr); 297421e10b8SBarry Smith low = i; 298421e10b8SBarry Smith } 299421e10b8SBarry Smith ailen[brow] = nrow; 300421e10b8SBarry Smith } 301421e10b8SBarry Smith A->same_nonzero = PETSC_FALSE; 302421e10b8SBarry Smith PetscFunctionReturn(0); 303421e10b8SBarry Smith } 304421e10b8SBarry Smith 305421e10b8SBarry Smith #undef __FUNCT__ 3065bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_BlockMat" 307112444f4SShri Abhyankar PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer) 308a5973382SShri Abhyankar { 309a5973382SShri Abhyankar PetscErrorCode ierr; 310a5973382SShri Abhyankar Mat tmpA; 311a5973382SShri Abhyankar PetscInt i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0; 312a5973382SShri Abhyankar const PetscInt *cols; 313a5973382SShri Abhyankar const PetscScalar *values; 314ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE,notdone; 315a5973382SShri Abhyankar Mat_SeqAIJ *a; 316a5973382SShri Abhyankar Mat_BlockMat *amat; 317a5973382SShri Abhyankar 318a5973382SShri Abhyankar PetscFunctionBegin; 319a5973382SShri Abhyankar ierr = MatCreate(PETSC_COMM_SELF,&tmpA);CHKERRQ(ierr); 320a5973382SShri Abhyankar ierr = MatSetType(tmpA,MATSEQAIJ);CHKERRQ(ierr); 321112444f4SShri Abhyankar ierr = MatLoad_SeqAIJ(tmpA,viewer);CHKERRQ(ierr); 322a5973382SShri Abhyankar 323a5973382SShri Abhyankar ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr); 324a5973382SShri Abhyankar ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr); 325a5973382SShri Abhyankar ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr); 326*acfcf0e5SJed Brown ierr = PetscOptionsBool("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 327a5973382SShri Abhyankar ierr = PetscOptionsEnd();CHKERRQ(ierr); 328a5973382SShri Abhyankar 329a5973382SShri Abhyankar /* Determine number of nonzero blocks for each block row */ 330a5973382SShri Abhyankar a = (Mat_SeqAIJ*) tmpA->data; 331a5973382SShri Abhyankar mbs = m/bs; 332a5973382SShri Abhyankar ierr = PetscMalloc3(mbs,PetscInt,&lens,bs,PetscInt*,&ii,bs,PetscInt,&ilens);CHKERRQ(ierr); 333a5973382SShri Abhyankar ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr); 334a5973382SShri Abhyankar 335a5973382SShri Abhyankar for (i=0; i<mbs; i++) { 336a5973382SShri Abhyankar for (j=0; j<bs; j++) { 337a5973382SShri Abhyankar ii[j] = a->j + a->i[i*bs + j]; 338a5973382SShri Abhyankar ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 339a5973382SShri Abhyankar } 340a5973382SShri Abhyankar 341a5973382SShri Abhyankar currentcol = -1; 342a5973382SShri Abhyankar notdone = PETSC_TRUE; 343a5973382SShri Abhyankar while (PETSC_TRUE) { 344a5973382SShri Abhyankar notdone = PETSC_FALSE; 345a5973382SShri Abhyankar nextcol = 1000000000; 346a5973382SShri Abhyankar for (j=0; j<bs; j++) { 347a5973382SShri Abhyankar while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { 348a5973382SShri Abhyankar ii[j]++; 349a5973382SShri Abhyankar ilens[j]--; 350a5973382SShri Abhyankar } 351a5973382SShri Abhyankar if (ilens[j] > 0) { 352a5973382SShri Abhyankar notdone = PETSC_TRUE; 353a5973382SShri Abhyankar nextcol = PetscMin(nextcol,ii[j][0]/bs); 354a5973382SShri Abhyankar } 355a5973382SShri Abhyankar } 356a5973382SShri Abhyankar if (!notdone) break; 357a5973382SShri Abhyankar if (!flg || (nextcol >= i)) lens[i]++; 358a5973382SShri Abhyankar currentcol = nextcol; 359a5973382SShri Abhyankar } 360a5973382SShri Abhyankar } 361a5973382SShri Abhyankar 362a5973382SShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 363a5973382SShri Abhyankar ierr = MatSetSizes(newmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 364a5973382SShri Abhyankar } 365a5973382SShri Abhyankar ierr = MatBlockMatSetPreallocation(newmat,bs,0,lens);CHKERRQ(ierr); 366a5973382SShri Abhyankar if (flg) { 367a5973382SShri Abhyankar ierr = MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 368a5973382SShri Abhyankar } 369a5973382SShri Abhyankar amat = (Mat_BlockMat*)(newmat)->data; 370a5973382SShri Abhyankar 371a5973382SShri Abhyankar /* preallocate the submatrices */ 372a5973382SShri Abhyankar ierr = PetscMalloc(bs*sizeof(PetscInt),&llens);CHKERRQ(ierr); 373a5973382SShri Abhyankar for (i=0; i<mbs; i++) { /* loops for block rows */ 374a5973382SShri Abhyankar for (j=0; j<bs; j++) { 375a5973382SShri Abhyankar ii[j] = a->j + a->i[i*bs + j]; 376a5973382SShri Abhyankar ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 377a5973382SShri Abhyankar } 378a5973382SShri Abhyankar 379a5973382SShri Abhyankar currentcol = 1000000000; 380a5973382SShri Abhyankar for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */ 381a5973382SShri Abhyankar if (ilens[j] > 0) { 382a5973382SShri Abhyankar currentcol = PetscMin(currentcol,ii[j][0]/bs); 383a5973382SShri Abhyankar } 384a5973382SShri Abhyankar } 385a5973382SShri Abhyankar 386a5973382SShri Abhyankar notdone = PETSC_TRUE; 387a5973382SShri Abhyankar while (PETSC_TRUE) { /* loops over blocks in block row */ 388a5973382SShri Abhyankar 389a5973382SShri Abhyankar notdone = PETSC_FALSE; 390a5973382SShri Abhyankar nextcol = 1000000000; 391a5973382SShri Abhyankar ierr = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr); 392a5973382SShri Abhyankar for (j=0; j<bs; j++) { /* loop over rows in block */ 393a5973382SShri Abhyankar while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */ 394a5973382SShri Abhyankar ii[j]++; 395a5973382SShri Abhyankar ilens[j]--; 396a5973382SShri Abhyankar llens[j]++; 397a5973382SShri Abhyankar } 398a5973382SShri Abhyankar if (ilens[j] > 0) { 399a5973382SShri Abhyankar notdone = PETSC_TRUE; 400a5973382SShri Abhyankar nextcol = PetscMin(nextcol,ii[j][0]/bs); 401a5973382SShri Abhyankar } 402a5973382SShri Abhyankar } 403a5973382SShri Abhyankar if (cnt >= amat->maxnz) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt); 404a5973382SShri Abhyankar if (!flg || currentcol >= i) { 405a5973382SShri Abhyankar amat->j[cnt] = currentcol; 406a5973382SShri Abhyankar ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr); 407a5973382SShri Abhyankar } 408a5973382SShri Abhyankar 409a5973382SShri Abhyankar if (!notdone) break; 410a5973382SShri Abhyankar currentcol = nextcol; 411a5973382SShri Abhyankar } 412a5973382SShri Abhyankar amat->ilen[i] = lens[i]; 413a5973382SShri Abhyankar } 414a5973382SShri Abhyankar CHKMEMQ; 415a5973382SShri Abhyankar 416a5973382SShri Abhyankar ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr); 417a5973382SShri Abhyankar ierr = PetscFree(llens);CHKERRQ(ierr); 418a5973382SShri Abhyankar 419a5973382SShri Abhyankar /* copy over the matrix, one row at a time */ 420a5973382SShri Abhyankar for (i=0; i<m; i++) { 421a5973382SShri Abhyankar ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 422a5973382SShri Abhyankar ierr = MatSetValues(newmat,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr); 423a5973382SShri Abhyankar ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 424a5973382SShri Abhyankar } 425a5973382SShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 426a5973382SShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 427a5973382SShri Abhyankar PetscFunctionReturn(0); 428a5973382SShri Abhyankar } 429a5973382SShri Abhyankar 430a5973382SShri Abhyankar #undef __FUNCT__ 431ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat" 432ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer) 433ccb205f8SBarry Smith { 43464075487SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 435ccb205f8SBarry Smith PetscErrorCode ierr; 436ccb205f8SBarry Smith const char *name; 437ccb205f8SBarry Smith PetscViewerFormat format; 438ccb205f8SBarry Smith 439ccb205f8SBarry Smith PetscFunctionBegin; 440ccb205f8SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 441ccb205f8SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 442ccb205f8SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 443ccb205f8SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr); 44464075487SBarry Smith if (A->symmetric) { 4458c1ad04dSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");CHKERRQ(ierr); 44664075487SBarry Smith } 447ccb205f8SBarry Smith } 448ccb205f8SBarry Smith PetscFunctionReturn(0); 449ccb205f8SBarry Smith } 450421e10b8SBarry Smith 451421e10b8SBarry Smith #undef __FUNCT__ 452421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat" 453421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat) 454421e10b8SBarry Smith { 455421e10b8SBarry Smith PetscErrorCode ierr; 456421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 457ccb205f8SBarry Smith PetscInt i; 458421e10b8SBarry Smith 459421e10b8SBarry Smith PetscFunctionBegin; 460421e10b8SBarry Smith if (bmat->right) { 461421e10b8SBarry Smith ierr = VecDestroy(bmat->right);CHKERRQ(ierr); 462421e10b8SBarry Smith } 463421e10b8SBarry Smith if (bmat->left) { 464421e10b8SBarry Smith ierr = VecDestroy(bmat->left);CHKERRQ(ierr); 465421e10b8SBarry Smith } 466290bbb0aSBarry Smith if (bmat->middle) { 467290bbb0aSBarry Smith ierr = VecDestroy(bmat->middle);CHKERRQ(ierr); 468290bbb0aSBarry Smith } 4696e63c7a1SBarry Smith if (bmat->workb) { 4706e63c7a1SBarry Smith ierr = VecDestroy(bmat->workb);CHKERRQ(ierr); 4716e63c7a1SBarry Smith } 472ccb205f8SBarry Smith if (bmat->diags) { 473d0f46423SBarry Smith for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) { 474ccb205f8SBarry Smith if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);} 475ccb205f8SBarry Smith } 476ccb205f8SBarry Smith } 477ccb205f8SBarry Smith if (bmat->a) { 478ccb205f8SBarry Smith for (i=0; i<bmat->nz; i++) { 479ccb205f8SBarry Smith if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);} 480ccb205f8SBarry Smith } 481ccb205f8SBarry Smith } 482ccb205f8SBarry Smith ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 483421e10b8SBarry Smith ierr = PetscFree(bmat);CHKERRQ(ierr); 484421e10b8SBarry Smith PetscFunctionReturn(0); 485421e10b8SBarry Smith } 486421e10b8SBarry Smith 487421e10b8SBarry Smith #undef __FUNCT__ 488421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat" 489421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 490421e10b8SBarry Smith { 491421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 492421e10b8SBarry Smith PetscErrorCode ierr; 493421e10b8SBarry Smith PetscScalar *xx,*yy; 494d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 495421e10b8SBarry Smith Mat *aa; 496421e10b8SBarry Smith 497421e10b8SBarry Smith PetscFunctionBegin; 498ccb205f8SBarry Smith CHKMEMQ; 499421e10b8SBarry Smith /* 500421e10b8SBarry Smith Standard CSR multiply except each entry is a Mat 501421e10b8SBarry Smith */ 502421e10b8SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 503ccb205f8SBarry Smith 504421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 505421e10b8SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 506421e10b8SBarry Smith aj = bmat->j; 507421e10b8SBarry Smith aa = bmat->a; 508421e10b8SBarry Smith ii = bmat->i; 509421e10b8SBarry Smith for (i=0; i<m; i++) { 510421e10b8SBarry Smith jrow = ii[i]; 511ccb205f8SBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 512421e10b8SBarry Smith n = ii[i+1] - jrow; 513421e10b8SBarry Smith for (j=0; j<n; j++) { 514421e10b8SBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 515ccb205f8SBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 516421e10b8SBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 517421e10b8SBarry Smith jrow++; 518421e10b8SBarry Smith } 519421e10b8SBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 520421e10b8SBarry Smith } 521421e10b8SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 522421e10b8SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 523ccb205f8SBarry Smith CHKMEMQ; 524421e10b8SBarry Smith PetscFunctionReturn(0); 525421e10b8SBarry Smith } 526421e10b8SBarry Smith 527421e10b8SBarry Smith #undef __FUNCT__ 528290bbb0aSBarry Smith #define __FUNCT__ "MatMult_BlockMat_Symmetric" 529290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y) 530290bbb0aSBarry Smith { 531290bbb0aSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 532290bbb0aSBarry Smith PetscErrorCode ierr; 533290bbb0aSBarry Smith PetscScalar *xx,*yy; 534d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 535290bbb0aSBarry Smith Mat *aa; 536290bbb0aSBarry Smith 537290bbb0aSBarry Smith PetscFunctionBegin; 538290bbb0aSBarry Smith CHKMEMQ; 539290bbb0aSBarry Smith /* 540290bbb0aSBarry Smith Standard CSR multiply except each entry is a Mat 541290bbb0aSBarry Smith */ 542290bbb0aSBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 543290bbb0aSBarry Smith 544290bbb0aSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 545290bbb0aSBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 546290bbb0aSBarry Smith aj = bmat->j; 547290bbb0aSBarry Smith aa = bmat->a; 548290bbb0aSBarry Smith ii = bmat->i; 549290bbb0aSBarry Smith for (i=0; i<m; i++) { 550290bbb0aSBarry Smith jrow = ii[i]; 551290bbb0aSBarry Smith n = ii[i+1] - jrow; 552290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 553290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr); 554290bbb0aSBarry Smith /* if we ALWAYS required a diagonal entry then could remove this if test */ 555290bbb0aSBarry Smith if (aj[jrow] == i) { 556290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 557290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 558290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 559290bbb0aSBarry Smith jrow++; 560290bbb0aSBarry Smith n--; 561290bbb0aSBarry Smith } 562290bbb0aSBarry Smith for (j=0; j<n; j++) { 563290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); /* upper triangular part */ 564290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 565290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 566290bbb0aSBarry Smith 567290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr); /* lower triangular part */ 568290bbb0aSBarry Smith ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr); 569290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 570290bbb0aSBarry Smith jrow++; 571290bbb0aSBarry Smith } 572290bbb0aSBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 573290bbb0aSBarry Smith ierr = VecResetArray(bmat->middle);CHKERRQ(ierr); 574290bbb0aSBarry Smith } 575290bbb0aSBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 576290bbb0aSBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 577290bbb0aSBarry Smith CHKMEMQ; 578290bbb0aSBarry Smith PetscFunctionReturn(0); 579290bbb0aSBarry Smith } 580290bbb0aSBarry Smith 581290bbb0aSBarry Smith #undef __FUNCT__ 582421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat" 583421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 584421e10b8SBarry Smith { 585421e10b8SBarry Smith PetscFunctionBegin; 586421e10b8SBarry Smith PetscFunctionReturn(0); 587421e10b8SBarry Smith } 588421e10b8SBarry Smith 589421e10b8SBarry Smith #undef __FUNCT__ 590421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat" 591421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 592421e10b8SBarry Smith { 593421e10b8SBarry Smith PetscFunctionBegin; 594421e10b8SBarry Smith PetscFunctionReturn(0); 595421e10b8SBarry Smith } 596421e10b8SBarry Smith 597421e10b8SBarry Smith #undef __FUNCT__ 598421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat" 599421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 600421e10b8SBarry Smith { 601421e10b8SBarry Smith PetscFunctionBegin; 602421e10b8SBarry Smith PetscFunctionReturn(0); 603421e10b8SBarry Smith } 604421e10b8SBarry Smith 605ccb205f8SBarry Smith /* 606ccb205f8SBarry Smith Adds diagonal pointers to sparse matrix structure. 607ccb205f8SBarry Smith */ 608ccb205f8SBarry Smith #undef __FUNCT__ 609ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat" 610ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) 611ccb205f8SBarry Smith { 612ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 613ccb205f8SBarry Smith PetscErrorCode ierr; 614d0f46423SBarry Smith PetscInt i,j,mbs = A->rmap->n/A->rmap->bs; 615ccb205f8SBarry Smith 616ccb205f8SBarry Smith PetscFunctionBegin; 617ccb205f8SBarry Smith if (!a->diag) { 618ccb205f8SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 619ccb205f8SBarry Smith } 620ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 621ccb205f8SBarry Smith a->diag[i] = a->i[i+1]; 622ccb205f8SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 623ccb205f8SBarry Smith if (a->j[j] == i) { 624ccb205f8SBarry Smith a->diag[i] = j; 625ccb205f8SBarry Smith break; 626ccb205f8SBarry Smith } 627ccb205f8SBarry Smith } 628ccb205f8SBarry Smith } 629ccb205f8SBarry Smith PetscFunctionReturn(0); 630ccb205f8SBarry Smith } 631ccb205f8SBarry Smith 632ccb205f8SBarry Smith #undef __FUNCT__ 633ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat" 6344aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 635ccb205f8SBarry Smith { 636ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 637ccb205f8SBarry Smith Mat_SeqAIJ *c; 638ccb205f8SBarry Smith PetscErrorCode ierr; 639ccb205f8SBarry Smith PetscInt i,k,first,step,lensi,nrows,ncols; 640d2a212eaSBarry Smith PetscInt *j_new,*i_new,*aj = a->j,*ailen = a->ilen; 6411620fd73SBarry Smith PetscScalar *a_new; 642ccb205f8SBarry Smith Mat C,*aa = a->a; 643ace3abfcSBarry Smith PetscBool stride,equal; 644ccb205f8SBarry Smith 645ccb205f8SBarry Smith PetscFunctionBegin; 646ccb205f8SBarry Smith ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr); 647e32f2f54SBarry Smith if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for idential column and row indices"); 648ccb205f8SBarry Smith ierr = ISStride(iscol,&stride);CHKERRQ(ierr); 649e32f2f54SBarry Smith if (!stride) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices"); 650ccb205f8SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 651e32f2f54SBarry Smith if (step != A->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block"); 652ccb205f8SBarry Smith 653ccb205f8SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 654ccb205f8SBarry Smith ncols = nrows; 655ccb205f8SBarry Smith 656ccb205f8SBarry Smith /* create submatrix */ 657ccb205f8SBarry Smith if (scall == MAT_REUSE_MATRIX) { 658ccb205f8SBarry Smith PetscInt n_cols,n_rows; 659ccb205f8SBarry Smith C = *B; 660ccb205f8SBarry Smith ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr); 661e32f2f54SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 662ccb205f8SBarry Smith ierr = MatZeroEntries(C);CHKERRQ(ierr); 663ccb205f8SBarry Smith } else { 6647adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 665ccb205f8SBarry Smith ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 6666e63c7a1SBarry Smith if (A->symmetric) { 6676e63c7a1SBarry Smith ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 6686e63c7a1SBarry Smith } else { 669ccb205f8SBarry Smith ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr); 6706e63c7a1SBarry Smith } 6716e63c7a1SBarry Smith ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr); 6726e63c7a1SBarry Smith ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr); 673ccb205f8SBarry Smith } 674ccb205f8SBarry Smith c = (Mat_SeqAIJ*)C->data; 675ccb205f8SBarry Smith 676ccb205f8SBarry Smith /* loop over rows inserting into submatrix */ 677ccb205f8SBarry Smith a_new = c->a; 678ccb205f8SBarry Smith j_new = c->j; 679ccb205f8SBarry Smith i_new = c->i; 680ccb205f8SBarry Smith 681ccb205f8SBarry Smith for (i=0; i<nrows; i++) { 682ccb205f8SBarry Smith lensi = ailen[i]; 683ccb205f8SBarry Smith for (k=0; k<lensi; k++) { 684ccb205f8SBarry Smith *j_new++ = *aj++; 6851620fd73SBarry Smith ierr = MatGetValue(*aa++,first,first,a_new++);CHKERRQ(ierr); 686ccb205f8SBarry Smith } 687ccb205f8SBarry Smith i_new[i+1] = i_new[i] + lensi; 688ccb205f8SBarry Smith c->ilen[i] = lensi; 689ccb205f8SBarry Smith } 690ccb205f8SBarry Smith 691ccb205f8SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 692ccb205f8SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 693ccb205f8SBarry Smith *B = C; 694ccb205f8SBarry Smith PetscFunctionReturn(0); 695ccb205f8SBarry Smith } 696ccb205f8SBarry Smith 697421e10b8SBarry Smith #undef __FUNCT__ 698421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat" 699421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 700421e10b8SBarry Smith { 701421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 702421e10b8SBarry Smith PetscErrorCode ierr; 703421e10b8SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 704ccb205f8SBarry Smith PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0; 705421e10b8SBarry Smith Mat *aa = a->a,*ap; 706421e10b8SBarry Smith 707421e10b8SBarry Smith PetscFunctionBegin; 708421e10b8SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 709421e10b8SBarry Smith 710421e10b8SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 711421e10b8SBarry Smith for (i=1; i<m; i++) { 712421e10b8SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 713421e10b8SBarry Smith fshift += imax[i-1] - ailen[i-1]; 714421e10b8SBarry Smith rmax = PetscMax(rmax,ailen[i]); 715421e10b8SBarry Smith if (fshift) { 716421e10b8SBarry Smith ip = aj + ai[i] ; 717421e10b8SBarry Smith ap = aa + ai[i] ; 718421e10b8SBarry Smith N = ailen[i]; 719421e10b8SBarry Smith for (j=0; j<N; j++) { 720421e10b8SBarry Smith ip[j-fshift] = ip[j]; 721421e10b8SBarry Smith ap[j-fshift] = ap[j]; 722421e10b8SBarry Smith } 723421e10b8SBarry Smith } 724421e10b8SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 725421e10b8SBarry Smith } 726421e10b8SBarry Smith if (m) { 727421e10b8SBarry Smith fshift += imax[m-1] - ailen[m-1]; 728421e10b8SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 729421e10b8SBarry Smith } 730421e10b8SBarry Smith /* reset ilen and imax for each row */ 731421e10b8SBarry Smith for (i=0; i<m; i++) { 732421e10b8SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 733421e10b8SBarry Smith } 734421e10b8SBarry Smith a->nz = ai[m]; 735ccb205f8SBarry Smith for (i=0; i<a->nz; i++) { 736ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG) 737e32f2f54SBarry Smith if (!aa[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz); 738ccb205f8SBarry Smith #endif 739ccb205f8SBarry Smith ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 740ccb205f8SBarry Smith ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 741ccb205f8SBarry Smith } 742ccb205f8SBarry Smith CHKMEMQ; 743d0f46423SBarry 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); 744421e10b8SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); 745421e10b8SBarry Smith ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); 7468e58a170SBarry Smith A->info.mallocs += a->reallocs; 747421e10b8SBarry Smith a->reallocs = 0; 748421e10b8SBarry Smith A->info.nz_unneeded = (double)fshift; 749421e10b8SBarry Smith a->rmax = rmax; 750421e10b8SBarry Smith 751421e10b8SBarry Smith A->same_nonzero = PETSC_TRUE; 752ccb205f8SBarry Smith ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr); 753421e10b8SBarry Smith PetscFunctionReturn(0); 754421e10b8SBarry Smith } 755421e10b8SBarry Smith 756290bbb0aSBarry Smith #undef __FUNCT__ 757290bbb0aSBarry Smith #define __FUNCT__ "MatSetOption_BlockMat" 758ace3abfcSBarry Smith PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscBool flg) 759290bbb0aSBarry Smith { 760290bbb0aSBarry Smith PetscFunctionBegin; 7614e0d8c25SBarry Smith if (opt == MAT_SYMMETRIC && flg) { 76241f059aeSBarry Smith A->ops->sor = MatSOR_BlockMat_Symmetric; 763290bbb0aSBarry Smith A->ops->mult = MatMult_BlockMat_Symmetric; 764290bbb0aSBarry Smith } else { 765290bbb0aSBarry Smith PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]); 766290bbb0aSBarry Smith } 767290bbb0aSBarry Smith PetscFunctionReturn(0); 768290bbb0aSBarry Smith } 769290bbb0aSBarry Smith 770290bbb0aSBarry Smith 771421e10b8SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 772421e10b8SBarry Smith 0, 773421e10b8SBarry Smith 0, 774421e10b8SBarry Smith MatMult_BlockMat, 775421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat, 776421e10b8SBarry Smith MatMultTranspose_BlockMat, 777421e10b8SBarry Smith MatMultTransposeAdd_BlockMat, 778421e10b8SBarry Smith 0, 779421e10b8SBarry Smith 0, 780421e10b8SBarry Smith 0, 781421e10b8SBarry Smith /*10*/ 0, 782421e10b8SBarry Smith 0, 783421e10b8SBarry Smith 0, 78441f059aeSBarry Smith MatSOR_BlockMat, 785421e10b8SBarry Smith 0, 786421e10b8SBarry Smith /*15*/ 0, 787421e10b8SBarry Smith 0, 788421e10b8SBarry Smith 0, 789421e10b8SBarry Smith 0, 790421e10b8SBarry Smith 0, 791421e10b8SBarry Smith /*20*/ 0, 792421e10b8SBarry Smith MatAssemblyEnd_BlockMat, 793290bbb0aSBarry Smith MatSetOption_BlockMat, 794421e10b8SBarry Smith 0, 795d519adbfSMatthew Knepley /*24*/ 0, 796421e10b8SBarry Smith 0, 797421e10b8SBarry Smith 0, 798421e10b8SBarry Smith 0, 799421e10b8SBarry Smith 0, 800d519adbfSMatthew Knepley /*29*/ 0, 801421e10b8SBarry Smith 0, 802421e10b8SBarry Smith 0, 803421e10b8SBarry Smith 0, 804421e10b8SBarry Smith 0, 805d519adbfSMatthew Knepley /*34*/ 0, 806421e10b8SBarry Smith 0, 807421e10b8SBarry Smith 0, 808421e10b8SBarry Smith 0, 809421e10b8SBarry Smith 0, 810d519adbfSMatthew Knepley /*39*/ 0, 811421e10b8SBarry Smith 0, 812421e10b8SBarry Smith 0, 813421e10b8SBarry Smith 0, 814421e10b8SBarry Smith 0, 815d519adbfSMatthew Knepley /*44*/ 0, 816421e10b8SBarry Smith 0, 817421e10b8SBarry Smith 0, 818421e10b8SBarry Smith 0, 819421e10b8SBarry Smith 0, 820d519adbfSMatthew Knepley /*49*/ 0, 821421e10b8SBarry Smith 0, 822421e10b8SBarry Smith 0, 823421e10b8SBarry Smith 0, 824421e10b8SBarry Smith 0, 825d519adbfSMatthew Knepley /*54*/ 0, 826421e10b8SBarry Smith 0, 827421e10b8SBarry Smith 0, 828421e10b8SBarry Smith 0, 829421e10b8SBarry Smith 0, 830d519adbfSMatthew Knepley /*59*/ MatGetSubMatrix_BlockMat, 831421e10b8SBarry Smith MatDestroy_BlockMat, 832ccb205f8SBarry Smith MatView_BlockMat, 833421e10b8SBarry Smith 0, 834421e10b8SBarry Smith 0, 835d519adbfSMatthew Knepley /*64*/ 0, 836421e10b8SBarry Smith 0, 837421e10b8SBarry Smith 0, 838421e10b8SBarry Smith 0, 839421e10b8SBarry Smith 0, 840d519adbfSMatthew Knepley /*69*/ 0, 841421e10b8SBarry Smith 0, 842421e10b8SBarry Smith 0, 843421e10b8SBarry Smith 0, 844421e10b8SBarry Smith 0, 845d519adbfSMatthew Knepley /*74*/ 0, 846421e10b8SBarry Smith 0, 847421e10b8SBarry Smith 0, 848421e10b8SBarry Smith 0, 849421e10b8SBarry Smith 0, 850d519adbfSMatthew Knepley /*79*/ 0, 851421e10b8SBarry Smith 0, 852421e10b8SBarry Smith 0, 853421e10b8SBarry Smith 0, 8545bba2384SShri Abhyankar MatLoad_BlockMat, 855d519adbfSMatthew Knepley /*84*/ 0, 856421e10b8SBarry Smith 0, 857421e10b8SBarry Smith 0, 858421e10b8SBarry Smith 0, 859421e10b8SBarry Smith 0, 860d519adbfSMatthew Knepley /*89*/ 0, 861421e10b8SBarry Smith 0, 862421e10b8SBarry Smith 0, 863421e10b8SBarry Smith 0, 864421e10b8SBarry Smith 0, 865d519adbfSMatthew Knepley /*94*/ 0, 866421e10b8SBarry Smith 0, 867421e10b8SBarry Smith 0, 868a5973382SShri Abhyankar 0, 869a5973382SShri Abhyankar 0, 870a5973382SShri Abhyankar /*99*/ 0, 871a5973382SShri Abhyankar 0, 872a5973382SShri Abhyankar 0, 873a5973382SShri Abhyankar 0, 874a5973382SShri Abhyankar 0, 875a5973382SShri Abhyankar /*104*/0, 876a5973382SShri Abhyankar 0, 877a5973382SShri Abhyankar 0, 878a5973382SShri Abhyankar 0, 879a5973382SShri Abhyankar 0, 880a5973382SShri Abhyankar /*109*/0, 881a5973382SShri Abhyankar 0, 882a5973382SShri Abhyankar 0, 883a5973382SShri Abhyankar 0, 884a5973382SShri Abhyankar 0, 885a5973382SShri Abhyankar /*114*/0, 886a5973382SShri Abhyankar 0, 887a5973382SShri Abhyankar 0, 888a5973382SShri Abhyankar 0, 889a5973382SShri Abhyankar 0, 890a5973382SShri Abhyankar /*119*/0, 891a5973382SShri Abhyankar 0, 892a5973382SShri Abhyankar 0, 8935bba2384SShri Abhyankar 0 894a5973382SShri Abhyankar }; 895421e10b8SBarry Smith 8968cdf2d9bSBarry Smith #undef __FUNCT__ 8978cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation" 8988cdf2d9bSBarry Smith /*@C 8998cdf2d9bSBarry Smith MatBlockMatSetPreallocation - For good matrix assembly performance 9008cdf2d9bSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 9018cdf2d9bSBarry Smith (or the array nnz). By setting these parameters accurately, performance 9028cdf2d9bSBarry Smith during matrix assembly can be increased by more than a factor of 50. 9038cdf2d9bSBarry Smith 9048cdf2d9bSBarry Smith Collective on MPI_Comm 9058cdf2d9bSBarry Smith 9068cdf2d9bSBarry Smith Input Parameters: 9078cdf2d9bSBarry Smith + B - The matrix 9088cdf2d9bSBarry Smith . bs - size of each block in matrix 9098cdf2d9bSBarry Smith . nz - number of nonzeros per block row (same for all rows) 9108cdf2d9bSBarry Smith - nnz - array containing the number of nonzeros in the various block rows 9118cdf2d9bSBarry Smith (possibly different for each row) or PETSC_NULL 9128cdf2d9bSBarry Smith 9138cdf2d9bSBarry Smith Notes: 9148cdf2d9bSBarry Smith If nnz is given then nz is ignored 9158cdf2d9bSBarry Smith 9168cdf2d9bSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 9178cdf2d9bSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 9188cdf2d9bSBarry Smith allocation. For large problems you MUST preallocate memory or you 9198cdf2d9bSBarry Smith will get TERRIBLE performance, see the users' manual chapter on matrices. 9208cdf2d9bSBarry Smith 9218cdf2d9bSBarry Smith Level: intermediate 9228cdf2d9bSBarry Smith 9238cdf2d9bSBarry Smith .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues() 9248cdf2d9bSBarry Smith 9258cdf2d9bSBarry Smith @*/ 9268cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 9278cdf2d9bSBarry Smith { 9284ac538c5SBarry Smith PetscErrorCode ierr; 9298cdf2d9bSBarry Smith 9308cdf2d9bSBarry Smith PetscFunctionBegin; 9314ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatBlockMatSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));CHKERRQ(ierr); 9328cdf2d9bSBarry Smith PetscFunctionReturn(0); 9338cdf2d9bSBarry Smith } 9348cdf2d9bSBarry Smith 9358cdf2d9bSBarry Smith EXTERN_C_BEGIN 9368cdf2d9bSBarry Smith #undef __FUNCT__ 9378cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation_BlockMat" 9388cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz) 9398cdf2d9bSBarry Smith { 9408cdf2d9bSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 9418cdf2d9bSBarry Smith PetscErrorCode ierr; 9428cdf2d9bSBarry Smith PetscInt i; 9438cdf2d9bSBarry Smith 9448cdf2d9bSBarry Smith PetscFunctionBegin; 94534ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 94634ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 94734ef9618SShri Abhyankar ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 94834ef9618SShri Abhyankar ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 94934ef9618SShri Abhyankar 95065e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Block size given %D must be great than zero",bs); 951e32f2f54SBarry Smith if (A->rmap->n % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap->n); 952e32f2f54SBarry Smith if (A->cmap->n % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap->n); 9538cdf2d9bSBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 954e32f2f54SBarry Smith if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 9558cdf2d9bSBarry Smith if (nnz) { 956d0f46423SBarry Smith for (i=0; i<A->rmap->n/bs; i++) { 957e32f2f54SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 958e32f2f54SBarry Smith if (nnz[i] > A->cmap->n/bs) SETERRQ3(PETSC_COMM_SELF,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); 9598cdf2d9bSBarry Smith } 9608cdf2d9bSBarry Smith } 961d0f46423SBarry Smith A->rmap->bs = A->cmap->bs = bs; 962d0f46423SBarry Smith bmat->mbs = A->rmap->n/bs; 9638cdf2d9bSBarry Smith 9648cdf2d9bSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr); 9658cdf2d9bSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->middle);CHKERRQ(ierr); 9668cdf2d9bSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr); 9678cdf2d9bSBarry Smith 9682ee49352SLisandro Dalcin if (!bmat->imax) { 969d0f46423SBarry Smith ierr = PetscMalloc2(A->rmap->n,PetscInt,&bmat->imax,A->rmap->n,PetscInt,&bmat->ilen);CHKERRQ(ierr); 970d0f46423SBarry Smith ierr = PetscLogObjectMemory(A,2*A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 9712ee49352SLisandro Dalcin } 9728cdf2d9bSBarry Smith if (nnz) { 9738cdf2d9bSBarry Smith nz = 0; 974d0f46423SBarry Smith for (i=0; i<A->rmap->n/A->rmap->bs; i++) { 9758cdf2d9bSBarry Smith bmat->imax[i] = nnz[i]; 9768cdf2d9bSBarry Smith nz += nnz[i]; 9778cdf2d9bSBarry Smith } 9788cdf2d9bSBarry Smith } else { 979e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation"); 9808cdf2d9bSBarry Smith } 9818cdf2d9bSBarry Smith 9828cdf2d9bSBarry Smith /* bmat->ilen will count nonzeros in each row so far. */ 9838cdf2d9bSBarry Smith for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;} 9848cdf2d9bSBarry Smith 9858cdf2d9bSBarry Smith /* allocate the matrix space */ 9862ee49352SLisandro Dalcin ierr = MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 987d0f46423SBarry Smith ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap->n+1,PetscInt,&bmat->i);CHKERRQ(ierr); 988d0f46423SBarry Smith ierr = PetscLogObjectMemory(A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 9898cdf2d9bSBarry Smith bmat->i[0] = 0; 9908cdf2d9bSBarry Smith for (i=1; i<bmat->mbs+1; i++) { 9918cdf2d9bSBarry Smith bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1]; 9928cdf2d9bSBarry Smith } 9938cdf2d9bSBarry Smith bmat->singlemalloc = PETSC_TRUE; 9948cdf2d9bSBarry Smith bmat->free_a = PETSC_TRUE; 9958cdf2d9bSBarry Smith bmat->free_ij = PETSC_TRUE; 9968cdf2d9bSBarry Smith 9978cdf2d9bSBarry Smith bmat->nz = 0; 9988cdf2d9bSBarry Smith bmat->maxnz = nz; 9998cdf2d9bSBarry Smith A->info.nz_unneeded = (double)bmat->maxnz; 10008cdf2d9bSBarry Smith 10018cdf2d9bSBarry Smith PetscFunctionReturn(0); 10028cdf2d9bSBarry Smith } 10038cdf2d9bSBarry Smith EXTERN_C_END 10048cdf2d9bSBarry Smith 1005421e10b8SBarry Smith /*MC 1006421e10b8SBarry Smith MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 1007421e10b8SBarry Smith consisting of (usually) sparse blocks. 1008421e10b8SBarry Smith 1009421e10b8SBarry Smith Level: advanced 1010421e10b8SBarry Smith 1011421e10b8SBarry Smith .seealso: MatCreateBlockMat() 1012421e10b8SBarry Smith 1013421e10b8SBarry Smith M*/ 1014421e10b8SBarry Smith 1015421e10b8SBarry Smith EXTERN_C_BEGIN 1016421e10b8SBarry Smith #undef __FUNCT__ 1017421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat" 1018421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A) 1019421e10b8SBarry Smith { 1020421e10b8SBarry Smith Mat_BlockMat *b; 1021421e10b8SBarry Smith PetscErrorCode ierr; 1022421e10b8SBarry Smith 1023421e10b8SBarry Smith PetscFunctionBegin; 102438f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_BlockMat,&b);CHKERRQ(ierr); 1025421e10b8SBarry Smith A->data = (void*)b; 102638f2d2fdSLisandro Dalcin ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1027421e10b8SBarry Smith 1028421e10b8SBarry Smith A->assembled = PETSC_TRUE; 1029421e10b8SBarry Smith A->preallocated = PETSC_FALSE; 1030421e10b8SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr); 1031290bbb0aSBarry Smith 10328cdf2d9bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatBlockMatSetPreallocation_C", 10338cdf2d9bSBarry Smith "MatBlockMatSetPreallocation_BlockMat", 10348cdf2d9bSBarry Smith MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr); 10358cdf2d9bSBarry Smith 1036421e10b8SBarry Smith PetscFunctionReturn(0); 1037421e10b8SBarry Smith } 1038421e10b8SBarry Smith EXTERN_C_END 1039421e10b8SBarry Smith 1040421e10b8SBarry Smith #undef __FUNCT__ 1041421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat" 1042421e10b8SBarry Smith /*@C 1043421e10b8SBarry Smith MatCreateBlockMat - Creates a new matrix based sparse Mat storage 1044421e10b8SBarry Smith 1045421e10b8SBarry Smith Collective on MPI_Comm 1046421e10b8SBarry Smith 1047421e10b8SBarry Smith Input Parameters: 1048421e10b8SBarry Smith + comm - MPI communicator 1049421e10b8SBarry Smith . m - number of rows 1050421e10b8SBarry Smith . n - number of columns 10518cdf2d9bSBarry Smith . bs - size of each submatrix 10528cdf2d9bSBarry Smith . nz - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known) 10538cdf2d9bSBarry Smith - nnz - expected number of nonzers per block row if known (use PETSC_NULL otherwise) 1054421e10b8SBarry Smith 1055421e10b8SBarry Smith 1056421e10b8SBarry Smith Output Parameter: 1057421e10b8SBarry Smith . A - the matrix 1058421e10b8SBarry Smith 1059421e10b8SBarry Smith Level: intermediate 1060421e10b8SBarry Smith 1061421e10b8SBarry Smith PETSc requires that matrices and vectors being used for certain 1062421e10b8SBarry Smith operations are partitioned accordingly. For example, when 1063421e10b8SBarry Smith creating a bmat matrix, A, that supports parallel matrix-vector 1064421e10b8SBarry Smith products using MatMult(A,x,y) the user should set the number 1065421e10b8SBarry Smith of local matrix rows to be the number of local elements of the 1066421e10b8SBarry Smith corresponding result vector, y. Note that this is information is 1067421e10b8SBarry Smith required for use of the matrix interface routines, even though 1068421e10b8SBarry Smith the bmat matrix may not actually be physically partitioned. 1069421e10b8SBarry Smith For example, 1070421e10b8SBarry Smith 1071421e10b8SBarry Smith .keywords: matrix, bmat, create 1072421e10b8SBarry Smith 1073421e10b8SBarry Smith .seealso: MATBLOCKMAT 1074421e10b8SBarry Smith @*/ 10758cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A) 1076421e10b8SBarry Smith { 1077421e10b8SBarry Smith PetscErrorCode ierr; 1078421e10b8SBarry Smith 1079421e10b8SBarry Smith PetscFunctionBegin; 1080421e10b8SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1081421e10b8SBarry Smith ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1082421e10b8SBarry Smith ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr); 10838cdf2d9bSBarry Smith ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 1084421e10b8SBarry Smith PetscFunctionReturn(0); 1085421e10b8SBarry Smith } 1086421e10b8SBarry Smith 1087421e10b8SBarry Smith 1088421e10b8SBarry Smith 1089