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 #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__ 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; 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); 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); 646e63c7a1SBarry Smith ierr = VecGetArray(a->workb,(PetscScalar**)&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); 1246e63c7a1SBarry Smith ierr = VecRestoreArray(a->workb,(PetscScalar**)&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; 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); 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); 166ccb205f8SBarry Smith ierr = VecGetArray(bb,(PetscScalar**)&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); 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; 240*71fd2e92SBarry 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) 246d0f46423SBarry Smith if (row >= A->rmap->N) SETERRQ2(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) 257d0f46423SBarry 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); 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; 281421e10b8SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 282421e10b8SBarry 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__ 306421e10b8SBarry Smith #define __FUNCT__ "MatLoad_BlockMat" 307a313700dSBarry Smith PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, const MatType type,Mat *A) 308421e10b8SBarry Smith { 309421e10b8SBarry Smith PetscErrorCode ierr; 310421e10b8SBarry Smith Mat tmpA; 311a0e1a203SBarry Smith PetscInt i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0; 312421e10b8SBarry Smith const PetscInt *cols; 313421e10b8SBarry Smith const PetscScalar *values; 31490d69ab7SBarry Smith PetscTruth flg = PETSC_FALSE,notdone; 3158cdf2d9bSBarry Smith Mat_SeqAIJ *a; 316a0e1a203SBarry Smith Mat_BlockMat *amat; 317421e10b8SBarry Smith 318421e10b8SBarry Smith PetscFunctionBegin; 319421e10b8SBarry Smith ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr); 320421e10b8SBarry Smith 321421e10b8SBarry Smith ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr); 32277925062SSatish Balay ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr); 323421e10b8SBarry Smith ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr); 32490d69ab7SBarry Smith ierr = PetscOptionsTruth("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 3258cdf2d9bSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 3268cdf2d9bSBarry Smith 327a0e1a203SBarry Smith /* Determine number of nonzero blocks for each block row */ 3288cdf2d9bSBarry Smith a = (Mat_SeqAIJ*) tmpA->data; 3298cdf2d9bSBarry Smith mbs = m/bs; 3301b453117SMatthew Knepley ierr = PetscMalloc3(mbs,PetscInt,&lens,bs,PetscInt*,&ii,bs,PetscInt,&ilens);CHKERRQ(ierr); 3318cdf2d9bSBarry Smith ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr); 3328cdf2d9bSBarry Smith 3338cdf2d9bSBarry Smith for (i=0; i<mbs; i++) { 3348cdf2d9bSBarry Smith for (j=0; j<bs; j++) { 3358cdf2d9bSBarry Smith ii[j] = a->j + a->i[i*bs + j]; 3368cdf2d9bSBarry Smith ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 3378cdf2d9bSBarry Smith } 338a0e1a203SBarry Smith 339a0e1a203SBarry Smith currentcol = -1; 3408cdf2d9bSBarry Smith notdone = PETSC_TRUE; 3418cdf2d9bSBarry Smith while (PETSC_TRUE) { 3428cdf2d9bSBarry Smith notdone = PETSC_FALSE; 3438cdf2d9bSBarry Smith nextcol = 1000000000; 3448cdf2d9bSBarry Smith for (j=0; j<bs; j++) { 345a0e1a203SBarry Smith while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { 3468cdf2d9bSBarry Smith ii[j]++; 3478cdf2d9bSBarry Smith ilens[j]--; 3488cdf2d9bSBarry Smith } 3498cdf2d9bSBarry Smith if (ilens[j] > 0) { 3508cdf2d9bSBarry Smith notdone = PETSC_TRUE; 3518cdf2d9bSBarry Smith nextcol = PetscMin(nextcol,ii[j][0]/bs); 3528cdf2d9bSBarry Smith } 3538cdf2d9bSBarry Smith } 3548cdf2d9bSBarry Smith if (!notdone) break; 355a0e1a203SBarry Smith if (!flg || (nextcol >= i)) lens[i]++; 3568cdf2d9bSBarry Smith currentcol = nextcol; 3578cdf2d9bSBarry Smith } 3588cdf2d9bSBarry Smith } 3598cdf2d9bSBarry Smith 3608cdf2d9bSBarry Smith ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,0,lens,A);CHKERRQ(ierr); 361290bbb0aSBarry Smith if (flg) { 3624e0d8c25SBarry Smith ierr = MatSetOption(*A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 363290bbb0aSBarry Smith } 364a0e1a203SBarry Smith amat = (Mat_BlockMat*)(*A)->data; 365a0e1a203SBarry Smith 366a0e1a203SBarry Smith /* preallocate the submatrices */ 367a0e1a203SBarry Smith ierr = PetscMalloc(bs*sizeof(PetscInt),&llens);CHKERRQ(ierr); 368a0e1a203SBarry Smith for (i=0; i<mbs; i++) { /* loops for block rows */ 369a0e1a203SBarry Smith for (j=0; j<bs; j++) { 370a0e1a203SBarry Smith ii[j] = a->j + a->i[i*bs + j]; 371a0e1a203SBarry Smith ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 372a0e1a203SBarry Smith } 373a0e1a203SBarry Smith 374a0e1a203SBarry Smith currentcol = 1000000000; 375a0e1a203SBarry Smith for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */ 376a0e1a203SBarry Smith if (ilens[j] > 0) { 377a0e1a203SBarry Smith currentcol = PetscMin(currentcol,ii[j][0]/bs); 378a0e1a203SBarry Smith } 379a0e1a203SBarry Smith } 380a0e1a203SBarry Smith 381a0e1a203SBarry Smith notdone = PETSC_TRUE; 382a0e1a203SBarry Smith while (PETSC_TRUE) { /* loops over blocks in block row */ 383a0e1a203SBarry Smith 384a0e1a203SBarry Smith notdone = PETSC_FALSE; 385a0e1a203SBarry Smith nextcol = 1000000000; 386a0e1a203SBarry Smith ierr = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr); 387a0e1a203SBarry Smith for (j=0; j<bs; j++) { /* loop over rows in block */ 388a0e1a203SBarry Smith while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */ 389a0e1a203SBarry Smith ii[j]++; 390a0e1a203SBarry Smith ilens[j]--; 391a0e1a203SBarry Smith llens[j]++; 392a0e1a203SBarry Smith } 393a0e1a203SBarry Smith if (ilens[j] > 0) { 394a0e1a203SBarry Smith notdone = PETSC_TRUE; 395a0e1a203SBarry Smith nextcol = PetscMin(nextcol,ii[j][0]/bs); 396a0e1a203SBarry Smith } 397a0e1a203SBarry Smith } 398a0e1a203SBarry Smith if (cnt >= amat->maxnz) SETERRQ1(PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt); 399a0e1a203SBarry Smith if (!flg || currentcol >= i) { 400a0e1a203SBarry Smith amat->j[cnt] = currentcol; 401a0e1a203SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr); 402a0e1a203SBarry Smith } 403a0e1a203SBarry Smith 404a0e1a203SBarry Smith if (!notdone) break; 405a0e1a203SBarry Smith currentcol = nextcol; 406a0e1a203SBarry Smith } 407a0e1a203SBarry Smith amat->ilen[i] = lens[i]; 408a0e1a203SBarry Smith } 409a0e1a203SBarry Smith CHKMEMQ; 410a0e1a203SBarry Smith 411a0e1a203SBarry Smith ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr); 412a0e1a203SBarry Smith ierr = PetscFree(llens);CHKERRQ(ierr); 413a0e1a203SBarry Smith 414a0e1a203SBarry Smith /* copy over the matrix, one row at a time */ 415421e10b8SBarry Smith for (i=0; i<m; i++) { 416421e10b8SBarry Smith ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 417421e10b8SBarry Smith ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr); 418ccb205f8SBarry Smith ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 419421e10b8SBarry Smith } 420421e10b8SBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 421421e10b8SBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 422421e10b8SBarry Smith PetscFunctionReturn(0); 423421e10b8SBarry Smith } 424421e10b8SBarry Smith 425ccb205f8SBarry Smith #undef __FUNCT__ 426ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat" 427ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer) 428ccb205f8SBarry Smith { 42964075487SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 430ccb205f8SBarry Smith PetscErrorCode ierr; 431ccb205f8SBarry Smith const char *name; 432ccb205f8SBarry Smith PetscViewerFormat format; 433ccb205f8SBarry Smith 434ccb205f8SBarry Smith PetscFunctionBegin; 435ccb205f8SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 436ccb205f8SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 437ccb205f8SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 438ccb205f8SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr); 43964075487SBarry Smith if (A->symmetric) { 4408c1ad04dSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");CHKERRQ(ierr); 44164075487SBarry Smith } 442ccb205f8SBarry Smith } 443ccb205f8SBarry Smith PetscFunctionReturn(0); 444ccb205f8SBarry Smith } 445421e10b8SBarry Smith 446421e10b8SBarry Smith #undef __FUNCT__ 447421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat" 448421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat) 449421e10b8SBarry Smith { 450421e10b8SBarry Smith PetscErrorCode ierr; 451421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 452ccb205f8SBarry Smith PetscInt i; 453421e10b8SBarry Smith 454421e10b8SBarry Smith PetscFunctionBegin; 455421e10b8SBarry Smith if (bmat->right) { 456421e10b8SBarry Smith ierr = VecDestroy(bmat->right);CHKERRQ(ierr); 457421e10b8SBarry Smith } 458421e10b8SBarry Smith if (bmat->left) { 459421e10b8SBarry Smith ierr = VecDestroy(bmat->left);CHKERRQ(ierr); 460421e10b8SBarry Smith } 461290bbb0aSBarry Smith if (bmat->middle) { 462290bbb0aSBarry Smith ierr = VecDestroy(bmat->middle);CHKERRQ(ierr); 463290bbb0aSBarry Smith } 4646e63c7a1SBarry Smith if (bmat->workb) { 4656e63c7a1SBarry Smith ierr = VecDestroy(bmat->workb);CHKERRQ(ierr); 4666e63c7a1SBarry Smith } 467ccb205f8SBarry Smith if (bmat->diags) { 468d0f46423SBarry Smith for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) { 469ccb205f8SBarry Smith if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);} 470ccb205f8SBarry Smith } 471ccb205f8SBarry Smith } 472ccb205f8SBarry Smith if (bmat->a) { 473ccb205f8SBarry Smith for (i=0; i<bmat->nz; i++) { 474ccb205f8SBarry Smith if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);} 475ccb205f8SBarry Smith } 476ccb205f8SBarry Smith } 477ccb205f8SBarry Smith ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 478421e10b8SBarry Smith ierr = PetscFree(bmat);CHKERRQ(ierr); 479421e10b8SBarry Smith PetscFunctionReturn(0); 480421e10b8SBarry Smith } 481421e10b8SBarry Smith 482421e10b8SBarry Smith #undef __FUNCT__ 483421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat" 484421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 485421e10b8SBarry Smith { 486421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 487421e10b8SBarry Smith PetscErrorCode ierr; 488421e10b8SBarry Smith PetscScalar *xx,*yy; 489d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 490421e10b8SBarry Smith Mat *aa; 491421e10b8SBarry Smith 492421e10b8SBarry Smith PetscFunctionBegin; 493ccb205f8SBarry Smith CHKMEMQ; 494421e10b8SBarry Smith /* 495421e10b8SBarry Smith Standard CSR multiply except each entry is a Mat 496421e10b8SBarry Smith */ 497421e10b8SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 498ccb205f8SBarry Smith 499421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 500421e10b8SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 501421e10b8SBarry Smith aj = bmat->j; 502421e10b8SBarry Smith aa = bmat->a; 503421e10b8SBarry Smith ii = bmat->i; 504421e10b8SBarry Smith for (i=0; i<m; i++) { 505421e10b8SBarry Smith jrow = ii[i]; 506ccb205f8SBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 507421e10b8SBarry Smith n = ii[i+1] - jrow; 508421e10b8SBarry Smith for (j=0; j<n; j++) { 509421e10b8SBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 510ccb205f8SBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 511421e10b8SBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 512421e10b8SBarry Smith jrow++; 513421e10b8SBarry Smith } 514421e10b8SBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 515421e10b8SBarry Smith } 516421e10b8SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 517421e10b8SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 518ccb205f8SBarry Smith CHKMEMQ; 519421e10b8SBarry Smith PetscFunctionReturn(0); 520421e10b8SBarry Smith } 521421e10b8SBarry Smith 522421e10b8SBarry Smith #undef __FUNCT__ 523290bbb0aSBarry Smith #define __FUNCT__ "MatMult_BlockMat_Symmetric" 524290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y) 525290bbb0aSBarry Smith { 526290bbb0aSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 527290bbb0aSBarry Smith PetscErrorCode ierr; 528290bbb0aSBarry Smith PetscScalar *xx,*yy; 529d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 530290bbb0aSBarry Smith Mat *aa; 531290bbb0aSBarry Smith 532290bbb0aSBarry Smith PetscFunctionBegin; 533290bbb0aSBarry Smith CHKMEMQ; 534290bbb0aSBarry Smith /* 535290bbb0aSBarry Smith Standard CSR multiply except each entry is a Mat 536290bbb0aSBarry Smith */ 537290bbb0aSBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 538290bbb0aSBarry Smith 539290bbb0aSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 540290bbb0aSBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 541290bbb0aSBarry Smith aj = bmat->j; 542290bbb0aSBarry Smith aa = bmat->a; 543290bbb0aSBarry Smith ii = bmat->i; 544290bbb0aSBarry Smith for (i=0; i<m; i++) { 545290bbb0aSBarry Smith jrow = ii[i]; 546290bbb0aSBarry Smith n = ii[i+1] - jrow; 547290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 548290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr); 549290bbb0aSBarry Smith /* if we ALWAYS required a diagonal entry then could remove this if test */ 550290bbb0aSBarry Smith if (aj[jrow] == i) { 551290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 552290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 553290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 554290bbb0aSBarry Smith jrow++; 555290bbb0aSBarry Smith n--; 556290bbb0aSBarry Smith } 557290bbb0aSBarry Smith for (j=0; j<n; j++) { 558290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); /* upper triangular part */ 559290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 560290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 561290bbb0aSBarry Smith 562290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr); /* lower triangular part */ 563290bbb0aSBarry Smith ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr); 564290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 565290bbb0aSBarry Smith jrow++; 566290bbb0aSBarry Smith } 567290bbb0aSBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 568290bbb0aSBarry Smith ierr = VecResetArray(bmat->middle);CHKERRQ(ierr); 569290bbb0aSBarry Smith } 570290bbb0aSBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 571290bbb0aSBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 572290bbb0aSBarry Smith CHKMEMQ; 573290bbb0aSBarry Smith PetscFunctionReturn(0); 574290bbb0aSBarry Smith } 575290bbb0aSBarry Smith 576290bbb0aSBarry Smith #undef __FUNCT__ 577421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat" 578421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 579421e10b8SBarry Smith { 580421e10b8SBarry Smith PetscFunctionBegin; 581421e10b8SBarry Smith PetscFunctionReturn(0); 582421e10b8SBarry Smith } 583421e10b8SBarry Smith 584421e10b8SBarry Smith #undef __FUNCT__ 585421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat" 586421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 587421e10b8SBarry Smith { 588421e10b8SBarry Smith PetscFunctionBegin; 589421e10b8SBarry Smith PetscFunctionReturn(0); 590421e10b8SBarry Smith } 591421e10b8SBarry Smith 592421e10b8SBarry Smith #undef __FUNCT__ 593421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat" 594421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 595421e10b8SBarry Smith { 596421e10b8SBarry Smith PetscFunctionBegin; 597421e10b8SBarry Smith PetscFunctionReturn(0); 598421e10b8SBarry Smith } 599421e10b8SBarry Smith 600ccb205f8SBarry Smith /* 601ccb205f8SBarry Smith Adds diagonal pointers to sparse matrix structure. 602ccb205f8SBarry Smith */ 603ccb205f8SBarry Smith #undef __FUNCT__ 604ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat" 605ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) 606ccb205f8SBarry Smith { 607ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 608ccb205f8SBarry Smith PetscErrorCode ierr; 609d0f46423SBarry Smith PetscInt i,j,mbs = A->rmap->n/A->rmap->bs; 610ccb205f8SBarry Smith 611ccb205f8SBarry Smith PetscFunctionBegin; 612ccb205f8SBarry Smith if (!a->diag) { 613ccb205f8SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 614ccb205f8SBarry Smith } 615ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 616ccb205f8SBarry Smith a->diag[i] = a->i[i+1]; 617ccb205f8SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 618ccb205f8SBarry Smith if (a->j[j] == i) { 619ccb205f8SBarry Smith a->diag[i] = j; 620ccb205f8SBarry Smith break; 621ccb205f8SBarry Smith } 622ccb205f8SBarry Smith } 623ccb205f8SBarry Smith } 624ccb205f8SBarry Smith PetscFunctionReturn(0); 625ccb205f8SBarry Smith } 626ccb205f8SBarry Smith 627ccb205f8SBarry Smith #undef __FUNCT__ 628ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat" 6294aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 630ccb205f8SBarry Smith { 631ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 632ccb205f8SBarry Smith Mat_SeqAIJ *c; 633ccb205f8SBarry Smith PetscErrorCode ierr; 634ccb205f8SBarry Smith PetscInt i,k,first,step,lensi,nrows,ncols; 635ccb205f8SBarry Smith PetscInt *j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen; 6361620fd73SBarry Smith PetscScalar *a_new; 637ccb205f8SBarry Smith Mat C,*aa = a->a; 638ccb205f8SBarry Smith PetscTruth stride,equal; 639ccb205f8SBarry Smith 640ccb205f8SBarry Smith PetscFunctionBegin; 641ccb205f8SBarry Smith ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr); 642ccb205f8SBarry Smith if (!equal) SETERRQ(PETSC_ERR_SUP,"Only for idential column and row indices"); 643ccb205f8SBarry Smith ierr = ISStride(iscol,&stride);CHKERRQ(ierr); 644ccb205f8SBarry Smith if (!stride) SETERRQ(PETSC_ERR_SUP,"Only for stride indices"); 645ccb205f8SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 646d0f46423SBarry Smith if (step != A->rmap->bs) SETERRQ(PETSC_ERR_SUP,"Can only select one entry from each block"); 647ccb205f8SBarry Smith 648ccb205f8SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 649ccb205f8SBarry Smith ncols = nrows; 650ccb205f8SBarry Smith 651ccb205f8SBarry Smith /* create submatrix */ 652ccb205f8SBarry Smith if (scall == MAT_REUSE_MATRIX) { 653ccb205f8SBarry Smith PetscInt n_cols,n_rows; 654ccb205f8SBarry Smith C = *B; 655ccb205f8SBarry Smith ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr); 656ccb205f8SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 657ccb205f8SBarry Smith ierr = MatZeroEntries(C);CHKERRQ(ierr); 658ccb205f8SBarry Smith } else { 6597adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 660ccb205f8SBarry Smith ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 6616e63c7a1SBarry Smith if (A->symmetric) { 6626e63c7a1SBarry Smith ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 6636e63c7a1SBarry Smith } else { 664ccb205f8SBarry Smith ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr); 6656e63c7a1SBarry Smith } 6666e63c7a1SBarry Smith ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr); 6676e63c7a1SBarry Smith ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr); 668ccb205f8SBarry Smith } 669ccb205f8SBarry Smith c = (Mat_SeqAIJ*)C->data; 670ccb205f8SBarry Smith 671ccb205f8SBarry Smith /* loop over rows inserting into submatrix */ 672ccb205f8SBarry Smith a_new = c->a; 673ccb205f8SBarry Smith j_new = c->j; 674ccb205f8SBarry Smith i_new = c->i; 675ccb205f8SBarry Smith 676ccb205f8SBarry Smith for (i=0; i<nrows; i++) { 677ccb205f8SBarry Smith ii = ai[i]; 678ccb205f8SBarry Smith lensi = ailen[i]; 679ccb205f8SBarry Smith for (k=0; k<lensi; k++) { 680ccb205f8SBarry Smith *j_new++ = *aj++; 6811620fd73SBarry Smith ierr = MatGetValue(*aa++,first,first,a_new++);CHKERRQ(ierr); 682ccb205f8SBarry Smith } 683ccb205f8SBarry Smith i_new[i+1] = i_new[i] + lensi; 684ccb205f8SBarry Smith c->ilen[i] = lensi; 685ccb205f8SBarry Smith } 686ccb205f8SBarry Smith 687ccb205f8SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 688ccb205f8SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 689ccb205f8SBarry Smith *B = C; 690ccb205f8SBarry Smith PetscFunctionReturn(0); 691ccb205f8SBarry Smith } 692ccb205f8SBarry Smith 693421e10b8SBarry Smith #undef __FUNCT__ 694421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat" 695421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 696421e10b8SBarry Smith { 697421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 698421e10b8SBarry Smith PetscErrorCode ierr; 699421e10b8SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 700ccb205f8SBarry Smith PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0; 701421e10b8SBarry Smith Mat *aa = a->a,*ap; 702421e10b8SBarry Smith 703421e10b8SBarry Smith PetscFunctionBegin; 704421e10b8SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 705421e10b8SBarry Smith 706421e10b8SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 707421e10b8SBarry Smith for (i=1; i<m; i++) { 708421e10b8SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 709421e10b8SBarry Smith fshift += imax[i-1] - ailen[i-1]; 710421e10b8SBarry Smith rmax = PetscMax(rmax,ailen[i]); 711421e10b8SBarry Smith if (fshift) { 712421e10b8SBarry Smith ip = aj + ai[i] ; 713421e10b8SBarry Smith ap = aa + ai[i] ; 714421e10b8SBarry Smith N = ailen[i]; 715421e10b8SBarry Smith for (j=0; j<N; j++) { 716421e10b8SBarry Smith ip[j-fshift] = ip[j]; 717421e10b8SBarry Smith ap[j-fshift] = ap[j]; 718421e10b8SBarry Smith } 719421e10b8SBarry Smith } 720421e10b8SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 721421e10b8SBarry Smith } 722421e10b8SBarry Smith if (m) { 723421e10b8SBarry Smith fshift += imax[m-1] - ailen[m-1]; 724421e10b8SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 725421e10b8SBarry Smith } 726421e10b8SBarry Smith /* reset ilen and imax for each row */ 727421e10b8SBarry Smith for (i=0; i<m; i++) { 728421e10b8SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 729421e10b8SBarry Smith } 730421e10b8SBarry Smith a->nz = ai[m]; 731ccb205f8SBarry Smith for (i=0; i<a->nz; i++) { 732ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG) 733ccb205f8SBarry Smith if (!aa[i]) SETERRQ3(PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz); 734ccb205f8SBarry Smith #endif 735ccb205f8SBarry Smith ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 736ccb205f8SBarry Smith ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 737ccb205f8SBarry Smith } 738ccb205f8SBarry Smith CHKMEMQ; 739d0f46423SBarry 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); 740421e10b8SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); 741421e10b8SBarry Smith ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); 742421e10b8SBarry Smith a->reallocs = 0; 743421e10b8SBarry Smith A->info.nz_unneeded = (double)fshift; 744421e10b8SBarry Smith a->rmax = rmax; 745421e10b8SBarry Smith 746421e10b8SBarry Smith A->same_nonzero = PETSC_TRUE; 747ccb205f8SBarry Smith ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr); 748421e10b8SBarry Smith PetscFunctionReturn(0); 749421e10b8SBarry Smith } 750421e10b8SBarry Smith 751290bbb0aSBarry Smith #undef __FUNCT__ 752290bbb0aSBarry Smith #define __FUNCT__ "MatSetOption_BlockMat" 7534e0d8c25SBarry Smith PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscTruth flg) 754290bbb0aSBarry Smith { 755290bbb0aSBarry Smith PetscFunctionBegin; 7564e0d8c25SBarry Smith if (opt == MAT_SYMMETRIC && flg) { 75741f059aeSBarry Smith A->ops->sor = MatSOR_BlockMat_Symmetric; 758290bbb0aSBarry Smith A->ops->mult = MatMult_BlockMat_Symmetric; 759290bbb0aSBarry Smith } else { 760290bbb0aSBarry Smith PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]); 761290bbb0aSBarry Smith } 762290bbb0aSBarry Smith PetscFunctionReturn(0); 763290bbb0aSBarry Smith } 764290bbb0aSBarry Smith 765290bbb0aSBarry Smith 766421e10b8SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 767421e10b8SBarry Smith 0, 768421e10b8SBarry Smith 0, 769421e10b8SBarry Smith MatMult_BlockMat, 770421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat, 771421e10b8SBarry Smith MatMultTranspose_BlockMat, 772421e10b8SBarry Smith MatMultTransposeAdd_BlockMat, 773421e10b8SBarry Smith 0, 774421e10b8SBarry Smith 0, 775421e10b8SBarry Smith 0, 776421e10b8SBarry Smith /*10*/ 0, 777421e10b8SBarry Smith 0, 778421e10b8SBarry Smith 0, 77941f059aeSBarry Smith MatSOR_BlockMat, 780421e10b8SBarry Smith 0, 781421e10b8SBarry Smith /*15*/ 0, 782421e10b8SBarry Smith 0, 783421e10b8SBarry Smith 0, 784421e10b8SBarry Smith 0, 785421e10b8SBarry Smith 0, 786421e10b8SBarry Smith /*20*/ 0, 787421e10b8SBarry Smith MatAssemblyEnd_BlockMat, 788290bbb0aSBarry Smith MatSetOption_BlockMat, 789421e10b8SBarry Smith 0, 790d519adbfSMatthew Knepley /*24*/ 0, 791421e10b8SBarry Smith 0, 792421e10b8SBarry Smith 0, 793421e10b8SBarry Smith 0, 794421e10b8SBarry Smith 0, 795d519adbfSMatthew Knepley /*29*/ 0, 796421e10b8SBarry Smith 0, 797421e10b8SBarry Smith 0, 798421e10b8SBarry Smith 0, 799421e10b8SBarry Smith 0, 800d519adbfSMatthew Knepley /*34*/ 0, 801421e10b8SBarry Smith 0, 802421e10b8SBarry Smith 0, 803421e10b8SBarry Smith 0, 804421e10b8SBarry Smith 0, 805d519adbfSMatthew Knepley /*39*/ 0, 806421e10b8SBarry Smith 0, 807421e10b8SBarry Smith 0, 808421e10b8SBarry Smith 0, 809421e10b8SBarry Smith 0, 810d519adbfSMatthew Knepley /*44*/ 0, 811421e10b8SBarry Smith 0, 812421e10b8SBarry Smith 0, 813421e10b8SBarry Smith 0, 814421e10b8SBarry Smith 0, 815d519adbfSMatthew Knepley /*49*/ 0, 816421e10b8SBarry Smith 0, 817421e10b8SBarry Smith 0, 818421e10b8SBarry Smith 0, 819421e10b8SBarry Smith 0, 820d519adbfSMatthew Knepley /*54*/ 0, 821421e10b8SBarry Smith 0, 822421e10b8SBarry Smith 0, 823421e10b8SBarry Smith 0, 824421e10b8SBarry Smith 0, 825d519adbfSMatthew Knepley /*59*/ MatGetSubMatrix_BlockMat, 826421e10b8SBarry Smith MatDestroy_BlockMat, 827ccb205f8SBarry Smith MatView_BlockMat, 828421e10b8SBarry Smith 0, 829421e10b8SBarry Smith 0, 830d519adbfSMatthew Knepley /*64*/ 0, 831421e10b8SBarry Smith 0, 832421e10b8SBarry Smith 0, 833421e10b8SBarry Smith 0, 834421e10b8SBarry Smith 0, 835d519adbfSMatthew Knepley /*69*/ 0, 836421e10b8SBarry Smith 0, 837421e10b8SBarry Smith 0, 838421e10b8SBarry Smith 0, 839421e10b8SBarry Smith 0, 840d519adbfSMatthew Knepley /*74*/ 0, 841421e10b8SBarry Smith 0, 842421e10b8SBarry Smith 0, 843421e10b8SBarry Smith 0, 844421e10b8SBarry Smith 0, 845d519adbfSMatthew Knepley /*79*/ 0, 846421e10b8SBarry Smith 0, 847421e10b8SBarry Smith 0, 848421e10b8SBarry Smith 0, 849421e10b8SBarry Smith MatLoad_BlockMat, 850d519adbfSMatthew Knepley /*84*/ 0, 851421e10b8SBarry Smith 0, 852421e10b8SBarry Smith 0, 853421e10b8SBarry Smith 0, 854421e10b8SBarry Smith 0, 855d519adbfSMatthew Knepley /*89*/ 0, 856421e10b8SBarry Smith 0, 857421e10b8SBarry Smith 0, 858421e10b8SBarry Smith 0, 859421e10b8SBarry Smith 0, 860d519adbfSMatthew Knepley /*94*/ 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 99526283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 99626283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 99726283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 99826283091SBarry Smith ierr = PetscLayoutSetUp(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