1421e10b8SBarry Smith #define PETSCMAT_DLL 2421e10b8SBarry Smith 3421e10b8SBarry Smith /* 4421e10b8SBarry Smith This provides a matrix that consists of Mats 5421e10b8SBarry Smith */ 6421e10b8SBarry Smith 7b9147fbbSdalcinl #include "include/private/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; 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; 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); 51*719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info);CHKERRQ(ierr); 52*719d5645SBarry 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); 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; 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; 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); 156*719d5645SBarry Smith ierr = MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info);CHKERRQ(ierr); 157*719d5645SBarry 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); 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; 232d0f46423SBarry 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) 245d0f46423SBarry 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) 256d0f46423SBarry 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 } 266421e10b8SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 267421e10b8SBarry Smith lastcol = col; 268421e10b8SBarry Smith while (high-low > 7) { 269421e10b8SBarry Smith t = (low+high)/2; 270421e10b8SBarry Smith if (rp[t] > bcol) high = t; 271421e10b8SBarry Smith else low = t; 272421e10b8SBarry Smith } 273421e10b8SBarry Smith for (i=low; i<high; i++) { 274421e10b8SBarry Smith if (rp[i] > bcol) break; 275421e10b8SBarry Smith if (rp[i] == bcol) { 276421e10b8SBarry Smith goto noinsert1; 277421e10b8SBarry Smith } 278421e10b8SBarry Smith } 279421e10b8SBarry Smith if (nonew == 1) goto noinsert1; 280421e10b8SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 281421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat); 282421e10b8SBarry Smith N = nrow++ - 1; high++; 283421e10b8SBarry Smith /* shift up all the later entries in this row */ 284421e10b8SBarry Smith for (ii=N; ii>=i; ii--) { 285421e10b8SBarry Smith rp[ii+1] = rp[ii]; 286421e10b8SBarry Smith ap[ii+1] = ap[ii]; 287421e10b8SBarry Smith } 288421e10b8SBarry Smith if (N>=i) ap[i] = 0; 289421e10b8SBarry Smith rp[i] = bcol; 290421e10b8SBarry Smith a->nz++; 291421e10b8SBarry Smith noinsert1:; 292421e10b8SBarry Smith if (!*(ap+i)) { 2936e63c7a1SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr); 294b0223207SBarry Smith } 295421e10b8SBarry Smith ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr); 296421e10b8SBarry Smith low = i; 297421e10b8SBarry Smith } 298421e10b8SBarry Smith ailen[brow] = nrow; 299421e10b8SBarry Smith } 300421e10b8SBarry Smith A->same_nonzero = PETSC_FALSE; 301421e10b8SBarry Smith PetscFunctionReturn(0); 302421e10b8SBarry Smith } 303421e10b8SBarry Smith 304421e10b8SBarry Smith #undef __FUNCT__ 305421e10b8SBarry Smith #define __FUNCT__ "MatLoad_BlockMat" 306a313700dSBarry Smith PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, const MatType type,Mat *A) 307421e10b8SBarry Smith { 308421e10b8SBarry Smith PetscErrorCode ierr; 309421e10b8SBarry Smith Mat tmpA; 310a0e1a203SBarry Smith PetscInt i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0; 311421e10b8SBarry Smith const PetscInt *cols; 312421e10b8SBarry Smith const PetscScalar *values; 3138cdf2d9bSBarry Smith PetscTruth flg,notdone; 3148cdf2d9bSBarry Smith Mat_SeqAIJ *a; 315a0e1a203SBarry Smith Mat_BlockMat *amat; 316421e10b8SBarry Smith 317421e10b8SBarry Smith PetscFunctionBegin; 318421e10b8SBarry Smith ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr); 319421e10b8SBarry Smith 320421e10b8SBarry Smith ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr); 32177925062SSatish Balay ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr); 322421e10b8SBarry Smith ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr); 323290bbb0aSBarry Smith ierr = PetscOptionsName("-matload_symmetric","Store the matrix as symmetric","MatLoad",&flg);CHKERRQ(ierr); 3248cdf2d9bSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 3258cdf2d9bSBarry Smith 326a0e1a203SBarry Smith /* Determine number of nonzero blocks for each block row */ 3278cdf2d9bSBarry Smith a = (Mat_SeqAIJ*) tmpA->data; 3288cdf2d9bSBarry Smith mbs = m/bs; 3291b453117SMatthew Knepley ierr = PetscMalloc3(mbs,PetscInt,&lens,bs,PetscInt*,&ii,bs,PetscInt,&ilens);CHKERRQ(ierr); 3308cdf2d9bSBarry Smith ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr); 3318cdf2d9bSBarry Smith 3328cdf2d9bSBarry Smith for (i=0; i<mbs; i++) { 3338cdf2d9bSBarry Smith for (j=0; j<bs; j++) { 3348cdf2d9bSBarry Smith ii[j] = a->j + a->i[i*bs + j]; 3358cdf2d9bSBarry Smith ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 3368cdf2d9bSBarry Smith } 337a0e1a203SBarry Smith 338a0e1a203SBarry Smith currentcol = -1; 3398cdf2d9bSBarry Smith notdone = PETSC_TRUE; 3408cdf2d9bSBarry Smith while (PETSC_TRUE) { 3418cdf2d9bSBarry Smith notdone = PETSC_FALSE; 3428cdf2d9bSBarry Smith nextcol = 1000000000; 3438cdf2d9bSBarry Smith for (j=0; j<bs; j++) { 344a0e1a203SBarry Smith while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { 3458cdf2d9bSBarry Smith ii[j]++; 3468cdf2d9bSBarry Smith ilens[j]--; 3478cdf2d9bSBarry Smith } 3488cdf2d9bSBarry Smith if (ilens[j] > 0) { 3498cdf2d9bSBarry Smith notdone = PETSC_TRUE; 3508cdf2d9bSBarry Smith nextcol = PetscMin(nextcol,ii[j][0]/bs); 3518cdf2d9bSBarry Smith } 3528cdf2d9bSBarry Smith } 3538cdf2d9bSBarry Smith if (!notdone) break; 354a0e1a203SBarry Smith if (!flg || (nextcol >= i)) lens[i]++; 3558cdf2d9bSBarry Smith currentcol = nextcol; 3568cdf2d9bSBarry Smith } 3578cdf2d9bSBarry Smith } 3588cdf2d9bSBarry Smith 3598cdf2d9bSBarry Smith ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,0,lens,A);CHKERRQ(ierr); 360290bbb0aSBarry Smith if (flg) { 3614e0d8c25SBarry Smith ierr = MatSetOption(*A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 362290bbb0aSBarry Smith } 363a0e1a203SBarry Smith amat = (Mat_BlockMat*)(*A)->data; 364a0e1a203SBarry Smith 365a0e1a203SBarry Smith /* preallocate the submatrices */ 366a0e1a203SBarry Smith ierr = PetscMalloc(bs*sizeof(PetscInt),&llens);CHKERRQ(ierr); 367a0e1a203SBarry Smith for (i=0; i<mbs; i++) { /* loops for block rows */ 368a0e1a203SBarry Smith for (j=0; j<bs; j++) { 369a0e1a203SBarry Smith ii[j] = a->j + a->i[i*bs + j]; 370a0e1a203SBarry Smith ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 371a0e1a203SBarry Smith } 372a0e1a203SBarry Smith 373a0e1a203SBarry Smith currentcol = 1000000000; 374a0e1a203SBarry Smith for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */ 375a0e1a203SBarry Smith if (ilens[j] > 0) { 376a0e1a203SBarry Smith currentcol = PetscMin(currentcol,ii[j][0]/bs); 377a0e1a203SBarry Smith } 378a0e1a203SBarry Smith } 379a0e1a203SBarry Smith 380a0e1a203SBarry Smith notdone = PETSC_TRUE; 381a0e1a203SBarry Smith while (PETSC_TRUE) { /* loops over blocks in block row */ 382a0e1a203SBarry Smith 383a0e1a203SBarry Smith notdone = PETSC_FALSE; 384a0e1a203SBarry Smith nextcol = 1000000000; 385a0e1a203SBarry Smith ierr = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr); 386a0e1a203SBarry Smith for (j=0; j<bs; j++) { /* loop over rows in block */ 387a0e1a203SBarry Smith while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */ 388a0e1a203SBarry Smith ii[j]++; 389a0e1a203SBarry Smith ilens[j]--; 390a0e1a203SBarry Smith llens[j]++; 391a0e1a203SBarry Smith } 392a0e1a203SBarry Smith if (ilens[j] > 0) { 393a0e1a203SBarry Smith notdone = PETSC_TRUE; 394a0e1a203SBarry Smith nextcol = PetscMin(nextcol,ii[j][0]/bs); 395a0e1a203SBarry Smith } 396a0e1a203SBarry Smith } 397a0e1a203SBarry Smith if (cnt >= amat->maxnz) SETERRQ1(PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt); 398a0e1a203SBarry Smith if (!flg || currentcol >= i) { 399a0e1a203SBarry Smith amat->j[cnt] = currentcol; 400a0e1a203SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr); 401a0e1a203SBarry Smith } 402a0e1a203SBarry Smith 403a0e1a203SBarry Smith if (!notdone) break; 404a0e1a203SBarry Smith currentcol = nextcol; 405a0e1a203SBarry Smith } 406a0e1a203SBarry Smith amat->ilen[i] = lens[i]; 407a0e1a203SBarry Smith } 408a0e1a203SBarry Smith CHKMEMQ; 409a0e1a203SBarry Smith 410a0e1a203SBarry Smith ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr); 411a0e1a203SBarry Smith ierr = PetscFree(llens);CHKERRQ(ierr); 412a0e1a203SBarry Smith 413a0e1a203SBarry Smith /* copy over the matrix, one row at a time */ 414421e10b8SBarry Smith for (i=0; i<m; i++) { 415421e10b8SBarry Smith ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 416421e10b8SBarry Smith ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr); 417ccb205f8SBarry Smith ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 418421e10b8SBarry Smith } 419421e10b8SBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 420421e10b8SBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 421421e10b8SBarry Smith PetscFunctionReturn(0); 422421e10b8SBarry Smith } 423421e10b8SBarry Smith 424ccb205f8SBarry Smith #undef __FUNCT__ 425ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat" 426ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer) 427ccb205f8SBarry Smith { 42864075487SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 429ccb205f8SBarry Smith PetscErrorCode ierr; 430ccb205f8SBarry Smith const char *name; 431ccb205f8SBarry Smith PetscViewerFormat format; 432ccb205f8SBarry Smith 433ccb205f8SBarry Smith PetscFunctionBegin; 434ccb205f8SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 435ccb205f8SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 436ccb205f8SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 437ccb205f8SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr); 43864075487SBarry Smith if (A->symmetric) { 4398c1ad04dSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");CHKERRQ(ierr); 44064075487SBarry Smith } 441ccb205f8SBarry Smith } 442ccb205f8SBarry Smith PetscFunctionReturn(0); 443ccb205f8SBarry Smith } 444421e10b8SBarry Smith 445421e10b8SBarry Smith #undef __FUNCT__ 446421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat" 447421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat) 448421e10b8SBarry Smith { 449421e10b8SBarry Smith PetscErrorCode ierr; 450421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 451ccb205f8SBarry Smith PetscInt i; 452421e10b8SBarry Smith 453421e10b8SBarry Smith PetscFunctionBegin; 454421e10b8SBarry Smith if (bmat->right) { 455421e10b8SBarry Smith ierr = VecDestroy(bmat->right);CHKERRQ(ierr); 456421e10b8SBarry Smith } 457421e10b8SBarry Smith if (bmat->left) { 458421e10b8SBarry Smith ierr = VecDestroy(bmat->left);CHKERRQ(ierr); 459421e10b8SBarry Smith } 460290bbb0aSBarry Smith if (bmat->middle) { 461290bbb0aSBarry Smith ierr = VecDestroy(bmat->middle);CHKERRQ(ierr); 462290bbb0aSBarry Smith } 4636e63c7a1SBarry Smith if (bmat->workb) { 4646e63c7a1SBarry Smith ierr = VecDestroy(bmat->workb);CHKERRQ(ierr); 4656e63c7a1SBarry Smith } 466ccb205f8SBarry Smith if (bmat->diags) { 467d0f46423SBarry Smith for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) { 468ccb205f8SBarry Smith if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);} 469ccb205f8SBarry Smith } 470ccb205f8SBarry Smith } 471ccb205f8SBarry Smith if (bmat->a) { 472ccb205f8SBarry Smith for (i=0; i<bmat->nz; i++) { 473ccb205f8SBarry Smith if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);} 474ccb205f8SBarry Smith } 475ccb205f8SBarry Smith } 476ccb205f8SBarry Smith ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 477421e10b8SBarry Smith ierr = PetscFree(bmat);CHKERRQ(ierr); 478421e10b8SBarry Smith PetscFunctionReturn(0); 479421e10b8SBarry Smith } 480421e10b8SBarry Smith 481421e10b8SBarry Smith #undef __FUNCT__ 482421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat" 483421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 484421e10b8SBarry Smith { 485421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 486421e10b8SBarry Smith PetscErrorCode ierr; 487421e10b8SBarry Smith PetscScalar *xx,*yy; 488d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 489421e10b8SBarry Smith Mat *aa; 490421e10b8SBarry Smith 491421e10b8SBarry Smith PetscFunctionBegin; 492ccb205f8SBarry Smith CHKMEMQ; 493421e10b8SBarry Smith /* 494421e10b8SBarry Smith Standard CSR multiply except each entry is a Mat 495421e10b8SBarry Smith */ 496421e10b8SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 497ccb205f8SBarry Smith 498421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 499421e10b8SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 500421e10b8SBarry Smith aj = bmat->j; 501421e10b8SBarry Smith aa = bmat->a; 502421e10b8SBarry Smith ii = bmat->i; 503421e10b8SBarry Smith for (i=0; i<m; i++) { 504421e10b8SBarry Smith jrow = ii[i]; 505ccb205f8SBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 506421e10b8SBarry Smith n = ii[i+1] - jrow; 507421e10b8SBarry Smith for (j=0; j<n; j++) { 508421e10b8SBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 509ccb205f8SBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 510421e10b8SBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 511421e10b8SBarry Smith jrow++; 512421e10b8SBarry Smith } 513421e10b8SBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 514421e10b8SBarry Smith } 515421e10b8SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 516421e10b8SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 517ccb205f8SBarry Smith CHKMEMQ; 518421e10b8SBarry Smith PetscFunctionReturn(0); 519421e10b8SBarry Smith } 520421e10b8SBarry Smith 521421e10b8SBarry Smith #undef __FUNCT__ 522290bbb0aSBarry Smith #define __FUNCT__ "MatMult_BlockMat_Symmetric" 523290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y) 524290bbb0aSBarry Smith { 525290bbb0aSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 526290bbb0aSBarry Smith PetscErrorCode ierr; 527290bbb0aSBarry Smith PetscScalar *xx,*yy; 528d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 529290bbb0aSBarry Smith Mat *aa; 530290bbb0aSBarry Smith 531290bbb0aSBarry Smith PetscFunctionBegin; 532290bbb0aSBarry Smith CHKMEMQ; 533290bbb0aSBarry Smith /* 534290bbb0aSBarry Smith Standard CSR multiply except each entry is a Mat 535290bbb0aSBarry Smith */ 536290bbb0aSBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 537290bbb0aSBarry Smith 538290bbb0aSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 539290bbb0aSBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 540290bbb0aSBarry Smith aj = bmat->j; 541290bbb0aSBarry Smith aa = bmat->a; 542290bbb0aSBarry Smith ii = bmat->i; 543290bbb0aSBarry Smith for (i=0; i<m; i++) { 544290bbb0aSBarry Smith jrow = ii[i]; 545290bbb0aSBarry Smith n = ii[i+1] - jrow; 546290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 547290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr); 548290bbb0aSBarry Smith /* if we ALWAYS required a diagonal entry then could remove this if test */ 549290bbb0aSBarry Smith if (aj[jrow] == i) { 550290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 551290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 552290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 553290bbb0aSBarry Smith jrow++; 554290bbb0aSBarry Smith n--; 555290bbb0aSBarry Smith } 556290bbb0aSBarry Smith for (j=0; j<n; j++) { 557290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); /* upper triangular part */ 558290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 559290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 560290bbb0aSBarry Smith 561290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr); /* lower triangular part */ 562290bbb0aSBarry Smith ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr); 563290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 564290bbb0aSBarry Smith jrow++; 565290bbb0aSBarry Smith } 566290bbb0aSBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 567290bbb0aSBarry Smith ierr = VecResetArray(bmat->middle);CHKERRQ(ierr); 568290bbb0aSBarry Smith } 569290bbb0aSBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 570290bbb0aSBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 571290bbb0aSBarry Smith CHKMEMQ; 572290bbb0aSBarry Smith PetscFunctionReturn(0); 573290bbb0aSBarry Smith } 574290bbb0aSBarry Smith 575290bbb0aSBarry Smith #undef __FUNCT__ 576421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat" 577421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 578421e10b8SBarry Smith { 579421e10b8SBarry Smith PetscFunctionBegin; 580421e10b8SBarry Smith PetscFunctionReturn(0); 581421e10b8SBarry Smith } 582421e10b8SBarry Smith 583421e10b8SBarry Smith #undef __FUNCT__ 584421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat" 585421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 586421e10b8SBarry Smith { 587421e10b8SBarry Smith PetscFunctionBegin; 588421e10b8SBarry Smith PetscFunctionReturn(0); 589421e10b8SBarry Smith } 590421e10b8SBarry Smith 591421e10b8SBarry Smith #undef __FUNCT__ 592421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat" 593421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 594421e10b8SBarry Smith { 595421e10b8SBarry Smith PetscFunctionBegin; 596421e10b8SBarry Smith PetscFunctionReturn(0); 597421e10b8SBarry Smith } 598421e10b8SBarry Smith 599ccb205f8SBarry Smith /* 600ccb205f8SBarry Smith Adds diagonal pointers to sparse matrix structure. 601ccb205f8SBarry Smith */ 602ccb205f8SBarry Smith #undef __FUNCT__ 603ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat" 604ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) 605ccb205f8SBarry Smith { 606ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 607ccb205f8SBarry Smith PetscErrorCode ierr; 608d0f46423SBarry Smith PetscInt i,j,mbs = A->rmap->n/A->rmap->bs; 609ccb205f8SBarry Smith 610ccb205f8SBarry Smith PetscFunctionBegin; 611ccb205f8SBarry Smith if (!a->diag) { 612ccb205f8SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 613ccb205f8SBarry Smith } 614ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 615ccb205f8SBarry Smith a->diag[i] = a->i[i+1]; 616ccb205f8SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 617ccb205f8SBarry Smith if (a->j[j] == i) { 618ccb205f8SBarry Smith a->diag[i] = j; 619ccb205f8SBarry Smith break; 620ccb205f8SBarry Smith } 621ccb205f8SBarry Smith } 622ccb205f8SBarry Smith } 623ccb205f8SBarry Smith PetscFunctionReturn(0); 624ccb205f8SBarry Smith } 625ccb205f8SBarry Smith 626ccb205f8SBarry Smith #undef __FUNCT__ 627ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat" 628ccb205f8SBarry Smith PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B) 629ccb205f8SBarry Smith { 630ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 631ccb205f8SBarry Smith Mat_SeqAIJ *c; 632ccb205f8SBarry Smith PetscErrorCode ierr; 633ccb205f8SBarry Smith PetscInt i,k,first,step,lensi,nrows,ncols; 634ccb205f8SBarry Smith PetscInt *j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen; 6351620fd73SBarry Smith PetscScalar *a_new; 636ccb205f8SBarry Smith Mat C,*aa = a->a; 637ccb205f8SBarry Smith PetscTruth stride,equal; 638ccb205f8SBarry Smith 639ccb205f8SBarry Smith PetscFunctionBegin; 640ccb205f8SBarry Smith ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr); 641ccb205f8SBarry Smith if (!equal) SETERRQ(PETSC_ERR_SUP,"Only for idential column and row indices"); 642ccb205f8SBarry Smith ierr = ISStride(iscol,&stride);CHKERRQ(ierr); 643ccb205f8SBarry Smith if (!stride) SETERRQ(PETSC_ERR_SUP,"Only for stride indices"); 644ccb205f8SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 645d0f46423SBarry Smith if (step != A->rmap->bs) SETERRQ(PETSC_ERR_SUP,"Can only select one entry from each block"); 646ccb205f8SBarry Smith 647ccb205f8SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 648ccb205f8SBarry Smith ncols = nrows; 649ccb205f8SBarry Smith 650ccb205f8SBarry Smith /* create submatrix */ 651ccb205f8SBarry Smith if (scall == MAT_REUSE_MATRIX) { 652ccb205f8SBarry Smith PetscInt n_cols,n_rows; 653ccb205f8SBarry Smith C = *B; 654ccb205f8SBarry Smith ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr); 655ccb205f8SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 656ccb205f8SBarry Smith ierr = MatZeroEntries(C);CHKERRQ(ierr); 657ccb205f8SBarry Smith } else { 6587adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 659ccb205f8SBarry Smith ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 6606e63c7a1SBarry Smith if (A->symmetric) { 6616e63c7a1SBarry Smith ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 6626e63c7a1SBarry Smith } else { 663ccb205f8SBarry Smith ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr); 6646e63c7a1SBarry Smith } 6656e63c7a1SBarry Smith ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr); 6666e63c7a1SBarry Smith ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr); 667ccb205f8SBarry Smith } 668ccb205f8SBarry Smith c = (Mat_SeqAIJ*)C->data; 669ccb205f8SBarry Smith 670ccb205f8SBarry Smith /* loop over rows inserting into submatrix */ 671ccb205f8SBarry Smith a_new = c->a; 672ccb205f8SBarry Smith j_new = c->j; 673ccb205f8SBarry Smith i_new = c->i; 674ccb205f8SBarry Smith 675ccb205f8SBarry Smith for (i=0; i<nrows; i++) { 676ccb205f8SBarry Smith ii = ai[i]; 677ccb205f8SBarry Smith lensi = ailen[i]; 678ccb205f8SBarry Smith for (k=0; k<lensi; k++) { 679ccb205f8SBarry Smith *j_new++ = *aj++; 6801620fd73SBarry Smith ierr = MatGetValue(*aa++,first,first,a_new++);CHKERRQ(ierr); 681ccb205f8SBarry Smith } 682ccb205f8SBarry Smith i_new[i+1] = i_new[i] + lensi; 683ccb205f8SBarry Smith c->ilen[i] = lensi; 684ccb205f8SBarry Smith } 685ccb205f8SBarry Smith 686ccb205f8SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 687ccb205f8SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 688ccb205f8SBarry Smith *B = C; 689ccb205f8SBarry Smith PetscFunctionReturn(0); 690ccb205f8SBarry Smith } 691ccb205f8SBarry Smith 692421e10b8SBarry Smith #undef __FUNCT__ 693421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat" 694421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 695421e10b8SBarry Smith { 696421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 697421e10b8SBarry Smith PetscErrorCode ierr; 698421e10b8SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 699ccb205f8SBarry Smith PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0; 700421e10b8SBarry Smith Mat *aa = a->a,*ap; 701421e10b8SBarry Smith 702421e10b8SBarry Smith PetscFunctionBegin; 703421e10b8SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 704421e10b8SBarry Smith 705421e10b8SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 706421e10b8SBarry Smith for (i=1; i<m; i++) { 707421e10b8SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 708421e10b8SBarry Smith fshift += imax[i-1] - ailen[i-1]; 709421e10b8SBarry Smith rmax = PetscMax(rmax,ailen[i]); 710421e10b8SBarry Smith if (fshift) { 711421e10b8SBarry Smith ip = aj + ai[i] ; 712421e10b8SBarry Smith ap = aa + ai[i] ; 713421e10b8SBarry Smith N = ailen[i]; 714421e10b8SBarry Smith for (j=0; j<N; j++) { 715421e10b8SBarry Smith ip[j-fshift] = ip[j]; 716421e10b8SBarry Smith ap[j-fshift] = ap[j]; 717421e10b8SBarry Smith } 718421e10b8SBarry Smith } 719421e10b8SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 720421e10b8SBarry Smith } 721421e10b8SBarry Smith if (m) { 722421e10b8SBarry Smith fshift += imax[m-1] - ailen[m-1]; 723421e10b8SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 724421e10b8SBarry Smith } 725421e10b8SBarry Smith /* reset ilen and imax for each row */ 726421e10b8SBarry Smith for (i=0; i<m; i++) { 727421e10b8SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 728421e10b8SBarry Smith } 729421e10b8SBarry Smith a->nz = ai[m]; 730ccb205f8SBarry Smith for (i=0; i<a->nz; i++) { 731ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG) 732ccb205f8SBarry Smith if (!aa[i]) SETERRQ3(PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz); 733ccb205f8SBarry Smith #endif 734ccb205f8SBarry Smith ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 735ccb205f8SBarry Smith ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 736ccb205f8SBarry Smith } 737ccb205f8SBarry Smith CHKMEMQ; 738d0f46423SBarry 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); 739421e10b8SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); 740421e10b8SBarry Smith ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); 741421e10b8SBarry Smith a->reallocs = 0; 742421e10b8SBarry Smith A->info.nz_unneeded = (double)fshift; 743421e10b8SBarry Smith a->rmax = rmax; 744421e10b8SBarry Smith 745421e10b8SBarry Smith A->same_nonzero = PETSC_TRUE; 746ccb205f8SBarry Smith ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr); 747421e10b8SBarry Smith PetscFunctionReturn(0); 748421e10b8SBarry Smith } 749421e10b8SBarry Smith 750290bbb0aSBarry Smith #undef __FUNCT__ 751290bbb0aSBarry Smith #define __FUNCT__ "MatSetOption_BlockMat" 7524e0d8c25SBarry Smith PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscTruth flg) 753290bbb0aSBarry Smith { 754290bbb0aSBarry Smith PetscFunctionBegin; 7554e0d8c25SBarry Smith if (opt == MAT_SYMMETRIC && flg) { 756290bbb0aSBarry Smith A->ops->relax = MatRelax_BlockMat_Symmetric; 757290bbb0aSBarry Smith A->ops->mult = MatMult_BlockMat_Symmetric; 758290bbb0aSBarry Smith } else { 759290bbb0aSBarry Smith PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]); 760290bbb0aSBarry Smith } 761290bbb0aSBarry Smith PetscFunctionReturn(0); 762290bbb0aSBarry Smith } 763290bbb0aSBarry Smith 764290bbb0aSBarry Smith 765421e10b8SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 766421e10b8SBarry Smith 0, 767421e10b8SBarry Smith 0, 768421e10b8SBarry Smith MatMult_BlockMat, 769421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat, 770421e10b8SBarry Smith MatMultTranspose_BlockMat, 771421e10b8SBarry Smith MatMultTransposeAdd_BlockMat, 772421e10b8SBarry Smith 0, 773421e10b8SBarry Smith 0, 774421e10b8SBarry Smith 0, 775421e10b8SBarry Smith /*10*/ 0, 776421e10b8SBarry Smith 0, 777421e10b8SBarry Smith 0, 778ccb205f8SBarry Smith MatRelax_BlockMat, 779421e10b8SBarry Smith 0, 780421e10b8SBarry Smith /*15*/ 0, 781421e10b8SBarry Smith 0, 782421e10b8SBarry Smith 0, 783421e10b8SBarry Smith 0, 784421e10b8SBarry Smith 0, 785421e10b8SBarry Smith /*20*/ 0, 786421e10b8SBarry Smith MatAssemblyEnd_BlockMat, 787421e10b8SBarry Smith 0, 788290bbb0aSBarry Smith MatSetOption_BlockMat, 789421e10b8SBarry Smith 0, 790421e10b8SBarry Smith /*25*/ 0, 791421e10b8SBarry Smith 0, 792421e10b8SBarry Smith 0, 793421e10b8SBarry Smith 0, 794421e10b8SBarry Smith 0, 795421e10b8SBarry Smith /*30*/ 0, 796421e10b8SBarry Smith 0, 797421e10b8SBarry Smith 0, 798421e10b8SBarry Smith 0, 799421e10b8SBarry Smith 0, 800421e10b8SBarry Smith /*35*/ 0, 801421e10b8SBarry Smith 0, 802421e10b8SBarry Smith 0, 803421e10b8SBarry Smith 0, 804421e10b8SBarry Smith 0, 805421e10b8SBarry Smith /*40*/ 0, 806421e10b8SBarry Smith 0, 807421e10b8SBarry Smith 0, 808421e10b8SBarry Smith 0, 809421e10b8SBarry Smith 0, 810421e10b8SBarry Smith /*45*/ 0, 811421e10b8SBarry Smith 0, 812421e10b8SBarry Smith 0, 813421e10b8SBarry Smith 0, 814421e10b8SBarry Smith 0, 8158cdf2d9bSBarry Smith /*50*/ 0, 816421e10b8SBarry Smith 0, 817421e10b8SBarry Smith 0, 818421e10b8SBarry Smith 0, 819421e10b8SBarry Smith 0, 820421e10b8SBarry Smith /*55*/ 0, 821421e10b8SBarry Smith 0, 822421e10b8SBarry Smith 0, 823421e10b8SBarry Smith 0, 824421e10b8SBarry Smith 0, 825ccb205f8SBarry Smith /*60*/ MatGetSubMatrix_BlockMat, 826421e10b8SBarry Smith MatDestroy_BlockMat, 827ccb205f8SBarry Smith MatView_BlockMat, 828421e10b8SBarry Smith 0, 829421e10b8SBarry Smith 0, 830421e10b8SBarry Smith /*65*/ 0, 831421e10b8SBarry Smith 0, 832421e10b8SBarry Smith 0, 833421e10b8SBarry Smith 0, 834421e10b8SBarry Smith 0, 835421e10b8SBarry Smith /*70*/ 0, 836421e10b8SBarry Smith 0, 837421e10b8SBarry Smith 0, 838421e10b8SBarry Smith 0, 839421e10b8SBarry Smith 0, 840421e10b8SBarry Smith /*75*/ 0, 841421e10b8SBarry Smith 0, 842421e10b8SBarry Smith 0, 843421e10b8SBarry Smith 0, 844421e10b8SBarry Smith 0, 845421e10b8SBarry Smith /*80*/ 0, 846421e10b8SBarry Smith 0, 847421e10b8SBarry Smith 0, 848421e10b8SBarry Smith 0, 849421e10b8SBarry Smith MatLoad_BlockMat, 850421e10b8SBarry Smith /*85*/ 0, 851421e10b8SBarry Smith 0, 852421e10b8SBarry Smith 0, 853421e10b8SBarry Smith 0, 854421e10b8SBarry Smith 0, 855421e10b8SBarry Smith /*90*/ 0, 856421e10b8SBarry Smith 0, 857421e10b8SBarry Smith 0, 858421e10b8SBarry Smith 0, 859421e10b8SBarry Smith 0, 860421e10b8SBarry Smith /*95*/ 0, 861421e10b8SBarry Smith 0, 862421e10b8SBarry Smith 0, 863421e10b8SBarry Smith 0}; 864421e10b8SBarry Smith 8658cdf2d9bSBarry Smith #undef __FUNCT__ 8668cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation" 8678cdf2d9bSBarry Smith /*@C 8688cdf2d9bSBarry Smith MatBlockMatSetPreallocation - For good matrix assembly performance 8698cdf2d9bSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 8708cdf2d9bSBarry Smith (or the array nnz). By setting these parameters accurately, performance 8718cdf2d9bSBarry Smith during matrix assembly can be increased by more than a factor of 50. 8728cdf2d9bSBarry Smith 8738cdf2d9bSBarry Smith Collective on MPI_Comm 8748cdf2d9bSBarry Smith 8758cdf2d9bSBarry Smith Input Parameters: 8768cdf2d9bSBarry Smith + B - The matrix 8778cdf2d9bSBarry Smith . bs - size of each block in matrix 8788cdf2d9bSBarry Smith . nz - number of nonzeros per block row (same for all rows) 8798cdf2d9bSBarry Smith - nnz - array containing the number of nonzeros in the various block rows 8808cdf2d9bSBarry Smith (possibly different for each row) or PETSC_NULL 8818cdf2d9bSBarry Smith 8828cdf2d9bSBarry Smith Notes: 8838cdf2d9bSBarry Smith If nnz is given then nz is ignored 8848cdf2d9bSBarry Smith 8858cdf2d9bSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 8868cdf2d9bSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 8878cdf2d9bSBarry Smith allocation. For large problems you MUST preallocate memory or you 8888cdf2d9bSBarry Smith will get TERRIBLE performance, see the users' manual chapter on matrices. 8898cdf2d9bSBarry Smith 8908cdf2d9bSBarry Smith Level: intermediate 8918cdf2d9bSBarry Smith 8928cdf2d9bSBarry Smith .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues() 8938cdf2d9bSBarry Smith 8948cdf2d9bSBarry Smith @*/ 8958cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 8968cdf2d9bSBarry Smith { 8978cdf2d9bSBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 8988cdf2d9bSBarry Smith 8998cdf2d9bSBarry Smith PetscFunctionBegin; 9008cdf2d9bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)B,"MatBlockMatSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 9018cdf2d9bSBarry Smith if (f) { 9028cdf2d9bSBarry Smith ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 9038cdf2d9bSBarry Smith } 9048cdf2d9bSBarry Smith PetscFunctionReturn(0); 9058cdf2d9bSBarry Smith } 9068cdf2d9bSBarry Smith 9078cdf2d9bSBarry Smith EXTERN_C_BEGIN 9088cdf2d9bSBarry Smith #undef __FUNCT__ 9098cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation_BlockMat" 9108cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz) 9118cdf2d9bSBarry Smith { 9128cdf2d9bSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 9138cdf2d9bSBarry Smith PetscErrorCode ierr; 9148cdf2d9bSBarry Smith PetscInt i; 9158cdf2d9bSBarry Smith 9168cdf2d9bSBarry Smith PetscFunctionBegin; 9178cdf2d9bSBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Block size given %D must be great than zero",bs); 918d0f46423SBarry Smith if (A->rmap->n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap->n); 919d0f46423SBarry Smith if (A->cmap->n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap->n); 9208cdf2d9bSBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 9218cdf2d9bSBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 9228cdf2d9bSBarry Smith if (nnz) { 923d0f46423SBarry Smith for (i=0; i<A->rmap->n/bs; i++) { 9248cdf2d9bSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 925d0f46423SBarry 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); 9268cdf2d9bSBarry Smith } 9278cdf2d9bSBarry Smith } 928d0f46423SBarry Smith A->rmap->bs = A->cmap->bs = bs; 929d0f46423SBarry Smith bmat->mbs = A->rmap->n/bs; 9308cdf2d9bSBarry Smith 9318cdf2d9bSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr); 9328cdf2d9bSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->middle);CHKERRQ(ierr); 9338cdf2d9bSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr); 9348cdf2d9bSBarry Smith 9352ee49352SLisandro Dalcin if (!bmat->imax) { 936d0f46423SBarry Smith ierr = PetscMalloc2(A->rmap->n,PetscInt,&bmat->imax,A->rmap->n,PetscInt,&bmat->ilen);CHKERRQ(ierr); 937d0f46423SBarry Smith ierr = PetscLogObjectMemory(A,2*A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 9382ee49352SLisandro Dalcin } 9398cdf2d9bSBarry Smith if (nnz) { 9408cdf2d9bSBarry Smith nz = 0; 941d0f46423SBarry Smith for (i=0; i<A->rmap->n/A->rmap->bs; i++) { 9428cdf2d9bSBarry Smith bmat->imax[i] = nnz[i]; 9438cdf2d9bSBarry Smith nz += nnz[i]; 9448cdf2d9bSBarry Smith } 9458cdf2d9bSBarry Smith } else { 9468cdf2d9bSBarry Smith SETERRQ(PETSC_ERR_SUP,"Currently requires block row by row preallocation"); 9478cdf2d9bSBarry Smith } 9488cdf2d9bSBarry Smith 9498cdf2d9bSBarry Smith /* bmat->ilen will count nonzeros in each row so far. */ 9508cdf2d9bSBarry Smith for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;} 9518cdf2d9bSBarry Smith 9528cdf2d9bSBarry Smith /* allocate the matrix space */ 9532ee49352SLisandro Dalcin ierr = MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 954d0f46423SBarry Smith ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap->n+1,PetscInt,&bmat->i);CHKERRQ(ierr); 955d0f46423SBarry Smith ierr = PetscLogObjectMemory(A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 9568cdf2d9bSBarry Smith bmat->i[0] = 0; 9578cdf2d9bSBarry Smith for (i=1; i<bmat->mbs+1; i++) { 9588cdf2d9bSBarry Smith bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1]; 9598cdf2d9bSBarry Smith } 9608cdf2d9bSBarry Smith bmat->singlemalloc = PETSC_TRUE; 9618cdf2d9bSBarry Smith bmat->free_a = PETSC_TRUE; 9628cdf2d9bSBarry Smith bmat->free_ij = PETSC_TRUE; 9638cdf2d9bSBarry Smith 9648cdf2d9bSBarry Smith bmat->nz = 0; 9658cdf2d9bSBarry Smith bmat->maxnz = nz; 9668cdf2d9bSBarry Smith A->info.nz_unneeded = (double)bmat->maxnz; 9678cdf2d9bSBarry Smith 9688cdf2d9bSBarry Smith PetscFunctionReturn(0); 9698cdf2d9bSBarry Smith } 9708cdf2d9bSBarry Smith EXTERN_C_END 9718cdf2d9bSBarry Smith 972421e10b8SBarry Smith /*MC 973421e10b8SBarry Smith MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 974421e10b8SBarry Smith consisting of (usually) sparse blocks. 975421e10b8SBarry Smith 976421e10b8SBarry Smith Level: advanced 977421e10b8SBarry Smith 978421e10b8SBarry Smith .seealso: MatCreateBlockMat() 979421e10b8SBarry Smith 980421e10b8SBarry Smith M*/ 981421e10b8SBarry Smith 982421e10b8SBarry Smith EXTERN_C_BEGIN 983421e10b8SBarry Smith #undef __FUNCT__ 984421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat" 985421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A) 986421e10b8SBarry Smith { 987421e10b8SBarry Smith Mat_BlockMat *b; 988421e10b8SBarry Smith PetscErrorCode ierr; 989421e10b8SBarry Smith 990421e10b8SBarry Smith PetscFunctionBegin; 99138f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_BlockMat,&b);CHKERRQ(ierr); 992421e10b8SBarry Smith A->data = (void*)b; 99338f2d2fdSLisandro Dalcin ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 994421e10b8SBarry Smith 995d0f46423SBarry Smith ierr = PetscMapSetBlockSize(A->rmap,1);CHKERRQ(ierr); 996d0f46423SBarry Smith ierr = PetscMapSetBlockSize(A->cmap,1);CHKERRQ(ierr); 997d0f46423SBarry Smith ierr = PetscMapSetUp(A->rmap);CHKERRQ(ierr); 998d0f46423SBarry Smith ierr = PetscMapSetUp(A->cmap);CHKERRQ(ierr); 999421e10b8SBarry Smith 1000421e10b8SBarry Smith A->assembled = PETSC_TRUE; 1001421e10b8SBarry Smith A->preallocated = PETSC_FALSE; 1002421e10b8SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr); 1003290bbb0aSBarry Smith 10048cdf2d9bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatBlockMatSetPreallocation_C", 10058cdf2d9bSBarry Smith "MatBlockMatSetPreallocation_BlockMat", 10068cdf2d9bSBarry Smith MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr); 10078cdf2d9bSBarry Smith 1008421e10b8SBarry Smith PetscFunctionReturn(0); 1009421e10b8SBarry Smith } 1010421e10b8SBarry Smith EXTERN_C_END 1011421e10b8SBarry Smith 1012421e10b8SBarry Smith #undef __FUNCT__ 1013421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat" 1014421e10b8SBarry Smith /*@C 1015421e10b8SBarry Smith MatCreateBlockMat - Creates a new matrix based sparse Mat storage 1016421e10b8SBarry Smith 1017421e10b8SBarry Smith Collective on MPI_Comm 1018421e10b8SBarry Smith 1019421e10b8SBarry Smith Input Parameters: 1020421e10b8SBarry Smith + comm - MPI communicator 1021421e10b8SBarry Smith . m - number of rows 1022421e10b8SBarry Smith . n - number of columns 10238cdf2d9bSBarry Smith . bs - size of each submatrix 10248cdf2d9bSBarry Smith . nz - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known) 10258cdf2d9bSBarry Smith - nnz - expected number of nonzers per block row if known (use PETSC_NULL otherwise) 1026421e10b8SBarry Smith 1027421e10b8SBarry Smith 1028421e10b8SBarry Smith Output Parameter: 1029421e10b8SBarry Smith . A - the matrix 1030421e10b8SBarry Smith 1031421e10b8SBarry Smith Level: intermediate 1032421e10b8SBarry Smith 1033421e10b8SBarry Smith PETSc requires that matrices and vectors being used for certain 1034421e10b8SBarry Smith operations are partitioned accordingly. For example, when 1035421e10b8SBarry Smith creating a bmat matrix, A, that supports parallel matrix-vector 1036421e10b8SBarry Smith products using MatMult(A,x,y) the user should set the number 1037421e10b8SBarry Smith of local matrix rows to be the number of local elements of the 1038421e10b8SBarry Smith corresponding result vector, y. Note that this is information is 1039421e10b8SBarry Smith required for use of the matrix interface routines, even though 1040421e10b8SBarry Smith the bmat matrix may not actually be physically partitioned. 1041421e10b8SBarry Smith For example, 1042421e10b8SBarry Smith 1043421e10b8SBarry Smith .keywords: matrix, bmat, create 1044421e10b8SBarry Smith 1045421e10b8SBarry Smith .seealso: MATBLOCKMAT 1046421e10b8SBarry Smith @*/ 10478cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A) 1048421e10b8SBarry Smith { 1049421e10b8SBarry Smith PetscErrorCode ierr; 1050421e10b8SBarry Smith 1051421e10b8SBarry Smith PetscFunctionBegin; 1052421e10b8SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1053421e10b8SBarry Smith ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1054421e10b8SBarry Smith ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr); 10558cdf2d9bSBarry Smith ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 1056421e10b8SBarry Smith PetscFunctionReturn(0); 1057421e10b8SBarry Smith } 1058421e10b8SBarry Smith 1059421e10b8SBarry Smith 1060421e10b8SBarry Smith 1061