1421e10b8SBarry Smith #define PETSCMAT_DLL 2421e10b8SBarry Smith 3421e10b8SBarry Smith /* 4421e10b8SBarry Smith This provides a matrix that consists of Mats 5421e10b8SBarry Smith */ 6421e10b8SBarry Smith 77c4f633dSBarry Smith #include "private/matimpl.h" /*I "petscmat.h" I*/ 87c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h" /* use the common AIJ data-structure */ 9ccb205f8SBarry Smith #include "petscksp.h" 10421e10b8SBarry Smith 11421e10b8SBarry Smith typedef struct { 12421e10b8SBarry Smith SEQAIJHEADER(Mat); 13421e10b8SBarry Smith SEQBAIJHEADER; 14ccb205f8SBarry Smith Mat *diags; 15421e10b8SBarry Smith 166e63c7a1SBarry Smith Vec left,right,middle,workb; /* dummy vectors to perform local parts of product */ 17421e10b8SBarry Smith } Mat_BlockMat; 18421e10b8SBarry Smith 19421e10b8SBarry Smith #undef __FUNCT__ 2041f059aeSBarry Smith #define __FUNCT__ "MatSOR_BlockMat_Symmetric" 2141f059aeSBarry Smith PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 22290bbb0aSBarry Smith { 23290bbb0aSBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 24290bbb0aSBarry Smith PetscScalar *x; 25290bbb0aSBarry Smith const Mat *v = a->a; 26290bbb0aSBarry Smith const PetscScalar *b; 27290bbb0aSBarry Smith PetscErrorCode ierr; 28d0f46423SBarry Smith PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs; 29290bbb0aSBarry Smith const PetscInt *idx; 30290bbb0aSBarry Smith IS row,col; 31290bbb0aSBarry Smith MatFactorInfo info; 326e63c7a1SBarry Smith Vec left = a->left,right = a->right, middle = a->middle; 33290bbb0aSBarry Smith Mat *diag; 34290bbb0aSBarry Smith 35290bbb0aSBarry Smith PetscFunctionBegin; 36290bbb0aSBarry Smith its = its*lits; 37e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 38e32f2f54SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 39e32f2f54SBarry Smith if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 40e32f2f54SBarry Smith if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift"); 41290bbb0aSBarry Smith if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) 42e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep"); 43290bbb0aSBarry Smith 44290bbb0aSBarry Smith if (!a->diags) { 45290bbb0aSBarry Smith ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr); 46290bbb0aSBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 47290bbb0aSBarry Smith for (i=0; i<mbs; i++) { 482692d6eeSBarry Smith ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr); 49719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info);CHKERRQ(ierr); 50719d5645SBarry Smith ierr = MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr); 51290bbb0aSBarry Smith ierr = ISDestroy(row);CHKERRQ(ierr); 52290bbb0aSBarry Smith ierr = ISDestroy(col);CHKERRQ(ierr); 53290bbb0aSBarry Smith } 546e63c7a1SBarry Smith ierr = VecDuplicate(bb,&a->workb);CHKERRQ(ierr); 55290bbb0aSBarry Smith } 56290bbb0aSBarry Smith diag = a->diags; 57290bbb0aSBarry Smith 58290bbb0aSBarry Smith ierr = VecSet(xx,0.0);CHKERRQ(ierr); 59290bbb0aSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 60290bbb0aSBarry Smith /* copy right hand side because it must be modified during iteration */ 616e63c7a1SBarry Smith ierr = VecCopy(bb,a->workb);CHKERRQ(ierr); 626e63c7a1SBarry Smith ierr = VecGetArray(a->workb,(PetscScalar**)&b);CHKERRQ(ierr); 63290bbb0aSBarry Smith 6441f059aeSBarry Smith /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */ 65290bbb0aSBarry Smith while (its--) { 66290bbb0aSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 67290bbb0aSBarry Smith 68290bbb0aSBarry Smith for (i=0; i<mbs; i++) { 696e63c7a1SBarry Smith n = a->i[i+1] - a->i[i] - 1; 706e63c7a1SBarry Smith idx = a->j + a->i[i] + 1; 716e63c7a1SBarry Smith v = a->a + a->i[i] + 1; 72290bbb0aSBarry Smith 73290bbb0aSBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 74290bbb0aSBarry Smith for (j=0; j<n; j++) { 75290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 76290bbb0aSBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 77290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 78290bbb0aSBarry Smith } 79290bbb0aSBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 80290bbb0aSBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 81290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 82290bbb0aSBarry Smith 83290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 84290bbb0aSBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 856e63c7a1SBarry Smith 8641f059aeSBarry Smith /* now adjust right hand side, see MatSOR_SeqSBAIJ */ 876e63c7a1SBarry Smith for (j=0; j<n; j++) { 886e63c7a1SBarry Smith ierr = MatMultTranspose(v[j],right,left);CHKERRQ(ierr); 896e63c7a1SBarry Smith ierr = VecPlaceArray(middle,b + idx[j]*bs);CHKERRQ(ierr); 906e63c7a1SBarry Smith ierr = VecAXPY(middle,-1.0,left);CHKERRQ(ierr); 916e63c7a1SBarry Smith ierr = VecResetArray(middle);CHKERRQ(ierr); 926e63c7a1SBarry Smith } 93290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 946e63c7a1SBarry Smith 95290bbb0aSBarry Smith } 96290bbb0aSBarry Smith } 97290bbb0aSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 98290bbb0aSBarry Smith 99290bbb0aSBarry Smith for (i=mbs-1; i>=0; i--) { 1006e63c7a1SBarry Smith n = a->i[i+1] - a->i[i] - 1; 1016e63c7a1SBarry Smith idx = a->j + a->i[i] + 1; 1026e63c7a1SBarry Smith v = a->a + a->i[i] + 1; 103290bbb0aSBarry Smith 104290bbb0aSBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 105290bbb0aSBarry Smith for (j=0; j<n; j++) { 106290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 107290bbb0aSBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 108290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 109290bbb0aSBarry Smith } 110290bbb0aSBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 111290bbb0aSBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 112290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 113290bbb0aSBarry Smith 114290bbb0aSBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 115290bbb0aSBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 116290bbb0aSBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 117290bbb0aSBarry Smith 118290bbb0aSBarry Smith } 119290bbb0aSBarry Smith } 120290bbb0aSBarry Smith } 121290bbb0aSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1226e63c7a1SBarry Smith ierr = VecRestoreArray(a->workb,(PetscScalar**)&b);CHKERRQ(ierr); 123290bbb0aSBarry Smith PetscFunctionReturn(0); 124290bbb0aSBarry Smith } 125290bbb0aSBarry Smith 126290bbb0aSBarry Smith #undef __FUNCT__ 12741f059aeSBarry Smith #define __FUNCT__ "MatSOR_BlockMat" 12841f059aeSBarry Smith PetscErrorCode MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 129ccb205f8SBarry Smith { 130ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 131ccb205f8SBarry Smith PetscScalar *x; 132ccb205f8SBarry Smith const Mat *v = a->a; 133ccb205f8SBarry Smith const PetscScalar *b; 134ccb205f8SBarry Smith PetscErrorCode ierr; 135d0f46423SBarry Smith PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs; 136ccb205f8SBarry Smith const PetscInt *idx; 137ccb205f8SBarry Smith IS row,col; 138ccb205f8SBarry Smith MatFactorInfo info; 139ccb205f8SBarry Smith Vec left = a->left,right = a->right; 140ccb205f8SBarry Smith Mat *diag; 141ccb205f8SBarry Smith 142ccb205f8SBarry Smith PetscFunctionBegin; 143ccb205f8SBarry Smith its = its*lits; 144e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 145e32f2f54SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 146e32f2f54SBarry Smith if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 147e32f2f54SBarry Smith if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift"); 148ccb205f8SBarry Smith 149ccb205f8SBarry Smith if (!a->diags) { 150ccb205f8SBarry Smith ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr); 151ccb205f8SBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 152ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 1532692d6eeSBarry Smith ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);CHKERRQ(ierr); 154719d5645SBarry Smith ierr = MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info);CHKERRQ(ierr); 155719d5645SBarry Smith ierr = MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);CHKERRQ(ierr); 156ccb205f8SBarry Smith ierr = ISDestroy(row);CHKERRQ(ierr); 157ccb205f8SBarry Smith ierr = ISDestroy(col);CHKERRQ(ierr); 158ccb205f8SBarry Smith } 159ccb205f8SBarry Smith } 160ccb205f8SBarry Smith diag = a->diags; 161ccb205f8SBarry Smith 162ccb205f8SBarry Smith ierr = VecSet(xx,0.0);CHKERRQ(ierr); 163ccb205f8SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 164ccb205f8SBarry Smith ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 165ccb205f8SBarry Smith 16641f059aeSBarry Smith /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */ 167ccb205f8SBarry Smith while (its--) { 168ccb205f8SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 169ccb205f8SBarry Smith 170ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 171ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 172ccb205f8SBarry Smith idx = a->j + a->i[i]; 173ccb205f8SBarry Smith v = a->a + a->i[i]; 174ccb205f8SBarry Smith 175ccb205f8SBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 176ccb205f8SBarry Smith for (j=0; j<n; j++) { 177ccb205f8SBarry Smith if (idx[j] != i) { 178ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 179ccb205f8SBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 180ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 181ccb205f8SBarry Smith } 182ccb205f8SBarry Smith } 183ccb205f8SBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 184ccb205f8SBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 185ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 186ccb205f8SBarry Smith 187ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 188ccb205f8SBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 189ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 190ccb205f8SBarry Smith } 191ccb205f8SBarry Smith } 192ccb205f8SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 193ccb205f8SBarry Smith 194ccb205f8SBarry Smith for (i=mbs-1; i>=0; i--) { 195ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 196ccb205f8SBarry Smith idx = a->j + a->i[i]; 197ccb205f8SBarry Smith v = a->a + a->i[i]; 198ccb205f8SBarry Smith 199ccb205f8SBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 200ccb205f8SBarry Smith for (j=0; j<n; j++) { 201ccb205f8SBarry Smith if (idx[j] != i) { 202ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 203ccb205f8SBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 204ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 205ccb205f8SBarry Smith } 206ccb205f8SBarry Smith } 207ccb205f8SBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 208ccb205f8SBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 209ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 210ccb205f8SBarry Smith 211ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 212ccb205f8SBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 213ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 214ccb205f8SBarry Smith 215ccb205f8SBarry Smith } 216ccb205f8SBarry Smith } 217ccb205f8SBarry Smith } 218ccb205f8SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 219ccb205f8SBarry Smith ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 220ccb205f8SBarry Smith PetscFunctionReturn(0); 221ccb205f8SBarry Smith } 222ccb205f8SBarry Smith 223ccb205f8SBarry Smith #undef __FUNCT__ 224421e10b8SBarry Smith #define __FUNCT__ "MatSetValues_BlockMat" 225421e10b8SBarry Smith PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 226421e10b8SBarry Smith { 227421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 228421e10b8SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 229421e10b8SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 230d0f46423SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol; 231421e10b8SBarry Smith PetscErrorCode ierr; 232421e10b8SBarry Smith PetscInt ridx,cidx; 233421e10b8SBarry Smith PetscTruth roworiented=a->roworiented; 234421e10b8SBarry Smith MatScalar value; 235421e10b8SBarry Smith Mat *ap,*aa = a->a; 236421e10b8SBarry Smith 237421e10b8SBarry Smith PetscFunctionBegin; 23871fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 239421e10b8SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 240421e10b8SBarry Smith row = im[k]; 241421e10b8SBarry Smith brow = row/bs; 242421e10b8SBarry Smith if (row < 0) continue; 243421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 244e32f2f54SBarry Smith if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1); 245421e10b8SBarry Smith #endif 246421e10b8SBarry Smith rp = aj + ai[brow]; 247421e10b8SBarry Smith ap = aa + ai[brow]; 248421e10b8SBarry Smith rmax = imax[brow]; 249421e10b8SBarry Smith nrow = ailen[brow]; 250421e10b8SBarry Smith low = 0; 251421e10b8SBarry Smith high = nrow; 252421e10b8SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 253421e10b8SBarry Smith if (in[l] < 0) continue; 254421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 255e32f2f54SBarry Smith if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1); 256421e10b8SBarry Smith #endif 257421e10b8SBarry Smith col = in[l]; bcol = col/bs; 2586e63c7a1SBarry Smith if (A->symmetric && brow > bcol) continue; 259421e10b8SBarry Smith ridx = row % bs; cidx = col % bs; 260421e10b8SBarry Smith if (roworiented) { 261421e10b8SBarry Smith value = v[l + k*n]; 262421e10b8SBarry Smith } else { 263421e10b8SBarry Smith value = v[k + l*m]; 264421e10b8SBarry Smith } 265421e10b8SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 266421e10b8SBarry Smith lastcol = col; 267421e10b8SBarry Smith while (high-low > 7) { 268421e10b8SBarry Smith t = (low+high)/2; 269421e10b8SBarry Smith if (rp[t] > bcol) high = t; 270421e10b8SBarry Smith else low = t; 271421e10b8SBarry Smith } 272421e10b8SBarry Smith for (i=low; i<high; i++) { 273421e10b8SBarry Smith if (rp[i] > bcol) break; 274421e10b8SBarry Smith if (rp[i] == bcol) { 275421e10b8SBarry Smith goto noinsert1; 276421e10b8SBarry Smith } 277421e10b8SBarry Smith } 278421e10b8SBarry Smith if (nonew == 1) goto noinsert1; 279e32f2f54SBarry Smith if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 280*fef13f97SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat); 281421e10b8SBarry Smith N = nrow++ - 1; high++; 282421e10b8SBarry Smith /* shift up all the later entries in this row */ 283421e10b8SBarry Smith for (ii=N; ii>=i; ii--) { 284421e10b8SBarry Smith rp[ii+1] = rp[ii]; 285421e10b8SBarry Smith ap[ii+1] = ap[ii]; 286421e10b8SBarry Smith } 287421e10b8SBarry Smith if (N>=i) ap[i] = 0; 288421e10b8SBarry Smith rp[i] = bcol; 289421e10b8SBarry Smith a->nz++; 290421e10b8SBarry Smith noinsert1:; 291421e10b8SBarry Smith if (!*(ap+i)) { 2926e63c7a1SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr); 293b0223207SBarry Smith } 294421e10b8SBarry Smith ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr); 295421e10b8SBarry Smith low = i; 296421e10b8SBarry Smith } 297421e10b8SBarry Smith ailen[brow] = nrow; 298421e10b8SBarry Smith } 299421e10b8SBarry Smith A->same_nonzero = PETSC_FALSE; 300421e10b8SBarry Smith PetscFunctionReturn(0); 301421e10b8SBarry Smith } 302421e10b8SBarry Smith 303421e10b8SBarry Smith #undef __FUNCT__ 304421e10b8SBarry Smith #define __FUNCT__ "MatLoad_BlockMat" 305a313700dSBarry Smith PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, const MatType type,Mat *A) 306421e10b8SBarry Smith { 307421e10b8SBarry Smith PetscErrorCode ierr; 308421e10b8SBarry Smith Mat tmpA; 309a0e1a203SBarry Smith PetscInt i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0; 310421e10b8SBarry Smith const PetscInt *cols; 311421e10b8SBarry Smith const PetscScalar *values; 31290d69ab7SBarry Smith PetscTruth flg = PETSC_FALSE,notdone; 3138cdf2d9bSBarry Smith Mat_SeqAIJ *a; 314a0e1a203SBarry Smith Mat_BlockMat *amat; 315421e10b8SBarry Smith 316421e10b8SBarry Smith PetscFunctionBegin; 317421e10b8SBarry Smith ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr); 318421e10b8SBarry Smith 319421e10b8SBarry Smith ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr); 32077925062SSatish Balay ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix 1","Mat");CHKERRQ(ierr); 321421e10b8SBarry Smith ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr); 32290d69ab7SBarry Smith ierr = PetscOptionsTruth("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 3238cdf2d9bSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 3248cdf2d9bSBarry Smith 325a0e1a203SBarry Smith /* Determine number of nonzero blocks for each block row */ 3268cdf2d9bSBarry Smith a = (Mat_SeqAIJ*) tmpA->data; 3278cdf2d9bSBarry Smith mbs = m/bs; 3281b453117SMatthew Knepley ierr = PetscMalloc3(mbs,PetscInt,&lens,bs,PetscInt*,&ii,bs,PetscInt,&ilens);CHKERRQ(ierr); 3298cdf2d9bSBarry Smith ierr = PetscMemzero(lens,mbs*sizeof(PetscInt));CHKERRQ(ierr); 3308cdf2d9bSBarry Smith 3318cdf2d9bSBarry Smith for (i=0; i<mbs; i++) { 3328cdf2d9bSBarry Smith for (j=0; j<bs; j++) { 3338cdf2d9bSBarry Smith ii[j] = a->j + a->i[i*bs + j]; 3348cdf2d9bSBarry Smith ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 3358cdf2d9bSBarry Smith } 336a0e1a203SBarry Smith 337a0e1a203SBarry Smith currentcol = -1; 3388cdf2d9bSBarry Smith notdone = PETSC_TRUE; 3398cdf2d9bSBarry Smith while (PETSC_TRUE) { 3408cdf2d9bSBarry Smith notdone = PETSC_FALSE; 3418cdf2d9bSBarry Smith nextcol = 1000000000; 3428cdf2d9bSBarry Smith for (j=0; j<bs; j++) { 343a0e1a203SBarry Smith while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { 3448cdf2d9bSBarry Smith ii[j]++; 3458cdf2d9bSBarry Smith ilens[j]--; 3468cdf2d9bSBarry Smith } 3478cdf2d9bSBarry Smith if (ilens[j] > 0) { 3488cdf2d9bSBarry Smith notdone = PETSC_TRUE; 3498cdf2d9bSBarry Smith nextcol = PetscMin(nextcol,ii[j][0]/bs); 3508cdf2d9bSBarry Smith } 3518cdf2d9bSBarry Smith } 3528cdf2d9bSBarry Smith if (!notdone) break; 353a0e1a203SBarry Smith if (!flg || (nextcol >= i)) lens[i]++; 3548cdf2d9bSBarry Smith currentcol = nextcol; 3558cdf2d9bSBarry Smith } 3568cdf2d9bSBarry Smith } 3578cdf2d9bSBarry Smith 3588cdf2d9bSBarry Smith ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,0,lens,A);CHKERRQ(ierr); 359290bbb0aSBarry Smith if (flg) { 3604e0d8c25SBarry Smith ierr = MatSetOption(*A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 361290bbb0aSBarry Smith } 362a0e1a203SBarry Smith amat = (Mat_BlockMat*)(*A)->data; 363a0e1a203SBarry Smith 364a0e1a203SBarry Smith /* preallocate the submatrices */ 365a0e1a203SBarry Smith ierr = PetscMalloc(bs*sizeof(PetscInt),&llens);CHKERRQ(ierr); 366a0e1a203SBarry Smith for (i=0; i<mbs; i++) { /* loops for block rows */ 367a0e1a203SBarry Smith for (j=0; j<bs; j++) { 368a0e1a203SBarry Smith ii[j] = a->j + a->i[i*bs + j]; 369a0e1a203SBarry Smith ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j]; 370a0e1a203SBarry Smith } 371a0e1a203SBarry Smith 372a0e1a203SBarry Smith currentcol = 1000000000; 373a0e1a203SBarry Smith for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */ 374a0e1a203SBarry Smith if (ilens[j] > 0) { 375a0e1a203SBarry Smith currentcol = PetscMin(currentcol,ii[j][0]/bs); 376a0e1a203SBarry Smith } 377a0e1a203SBarry Smith } 378a0e1a203SBarry Smith 379a0e1a203SBarry Smith notdone = PETSC_TRUE; 380a0e1a203SBarry Smith while (PETSC_TRUE) { /* loops over blocks in block row */ 381a0e1a203SBarry Smith 382a0e1a203SBarry Smith notdone = PETSC_FALSE; 383a0e1a203SBarry Smith nextcol = 1000000000; 384a0e1a203SBarry Smith ierr = PetscMemzero(llens,bs*sizeof(PetscInt));CHKERRQ(ierr); 385a0e1a203SBarry Smith for (j=0; j<bs; j++) { /* loop over rows in block */ 386a0e1a203SBarry Smith while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */ 387a0e1a203SBarry Smith ii[j]++; 388a0e1a203SBarry Smith ilens[j]--; 389a0e1a203SBarry Smith llens[j]++; 390a0e1a203SBarry Smith } 391a0e1a203SBarry Smith if (ilens[j] > 0) { 392a0e1a203SBarry Smith notdone = PETSC_TRUE; 393a0e1a203SBarry Smith nextcol = PetscMin(nextcol,ii[j][0]/bs); 394a0e1a203SBarry Smith } 395a0e1a203SBarry Smith } 396e32f2f54SBarry Smith if (cnt >= amat->maxnz) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt); 397a0e1a203SBarry Smith if (!flg || currentcol >= i) { 398a0e1a203SBarry Smith amat->j[cnt] = currentcol; 399a0e1a203SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);CHKERRQ(ierr); 400a0e1a203SBarry Smith } 401a0e1a203SBarry Smith 402a0e1a203SBarry Smith if (!notdone) break; 403a0e1a203SBarry Smith currentcol = nextcol; 404a0e1a203SBarry Smith } 405a0e1a203SBarry Smith amat->ilen[i] = lens[i]; 406a0e1a203SBarry Smith } 407a0e1a203SBarry Smith CHKMEMQ; 408a0e1a203SBarry Smith 409a0e1a203SBarry Smith ierr = PetscFree3(lens,ii,ilens);CHKERRQ(ierr); 410a0e1a203SBarry Smith ierr = PetscFree(llens);CHKERRQ(ierr); 411a0e1a203SBarry Smith 412a0e1a203SBarry Smith /* copy over the matrix, one row at a time */ 413421e10b8SBarry Smith for (i=0; i<m; i++) { 414421e10b8SBarry Smith ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 415421e10b8SBarry Smith ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr); 416ccb205f8SBarry Smith ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 417421e10b8SBarry Smith } 418421e10b8SBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 419421e10b8SBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 420421e10b8SBarry Smith PetscFunctionReturn(0); 421421e10b8SBarry Smith } 422421e10b8SBarry Smith 423ccb205f8SBarry Smith #undef __FUNCT__ 424ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat" 425ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer) 426ccb205f8SBarry Smith { 42764075487SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 428ccb205f8SBarry Smith PetscErrorCode ierr; 429ccb205f8SBarry Smith const char *name; 430ccb205f8SBarry Smith PetscViewerFormat format; 431ccb205f8SBarry Smith 432ccb205f8SBarry Smith PetscFunctionBegin; 433ccb205f8SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 434ccb205f8SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 435ccb205f8SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 436ccb205f8SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr); 43764075487SBarry Smith if (A->symmetric) { 4388c1ad04dSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");CHKERRQ(ierr); 43964075487SBarry Smith } 440ccb205f8SBarry Smith } 441ccb205f8SBarry Smith PetscFunctionReturn(0); 442ccb205f8SBarry Smith } 443421e10b8SBarry Smith 444421e10b8SBarry Smith #undef __FUNCT__ 445421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat" 446421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat) 447421e10b8SBarry Smith { 448421e10b8SBarry Smith PetscErrorCode ierr; 449421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 450ccb205f8SBarry Smith PetscInt i; 451421e10b8SBarry Smith 452421e10b8SBarry Smith PetscFunctionBegin; 453421e10b8SBarry Smith if (bmat->right) { 454421e10b8SBarry Smith ierr = VecDestroy(bmat->right);CHKERRQ(ierr); 455421e10b8SBarry Smith } 456421e10b8SBarry Smith if (bmat->left) { 457421e10b8SBarry Smith ierr = VecDestroy(bmat->left);CHKERRQ(ierr); 458421e10b8SBarry Smith } 459290bbb0aSBarry Smith if (bmat->middle) { 460290bbb0aSBarry Smith ierr = VecDestroy(bmat->middle);CHKERRQ(ierr); 461290bbb0aSBarry Smith } 4626e63c7a1SBarry Smith if (bmat->workb) { 4636e63c7a1SBarry Smith ierr = VecDestroy(bmat->workb);CHKERRQ(ierr); 4646e63c7a1SBarry Smith } 465ccb205f8SBarry Smith if (bmat->diags) { 466d0f46423SBarry Smith for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) { 467ccb205f8SBarry Smith if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);} 468ccb205f8SBarry Smith } 469ccb205f8SBarry Smith } 470ccb205f8SBarry Smith if (bmat->a) { 471ccb205f8SBarry Smith for (i=0; i<bmat->nz; i++) { 472ccb205f8SBarry Smith if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);} 473ccb205f8SBarry Smith } 474ccb205f8SBarry Smith } 475ccb205f8SBarry Smith ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 476421e10b8SBarry Smith ierr = PetscFree(bmat);CHKERRQ(ierr); 477421e10b8SBarry Smith PetscFunctionReturn(0); 478421e10b8SBarry Smith } 479421e10b8SBarry Smith 480421e10b8SBarry Smith #undef __FUNCT__ 481421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat" 482421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 483421e10b8SBarry Smith { 484421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 485421e10b8SBarry Smith PetscErrorCode ierr; 486421e10b8SBarry Smith PetscScalar *xx,*yy; 487d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 488421e10b8SBarry Smith Mat *aa; 489421e10b8SBarry Smith 490421e10b8SBarry Smith PetscFunctionBegin; 491ccb205f8SBarry Smith CHKMEMQ; 492421e10b8SBarry Smith /* 493421e10b8SBarry Smith Standard CSR multiply except each entry is a Mat 494421e10b8SBarry Smith */ 495421e10b8SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 496ccb205f8SBarry Smith 497421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 498421e10b8SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 499421e10b8SBarry Smith aj = bmat->j; 500421e10b8SBarry Smith aa = bmat->a; 501421e10b8SBarry Smith ii = bmat->i; 502421e10b8SBarry Smith for (i=0; i<m; i++) { 503421e10b8SBarry Smith jrow = ii[i]; 504ccb205f8SBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 505421e10b8SBarry Smith n = ii[i+1] - jrow; 506421e10b8SBarry Smith for (j=0; j<n; j++) { 507421e10b8SBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 508ccb205f8SBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 509421e10b8SBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 510421e10b8SBarry Smith jrow++; 511421e10b8SBarry Smith } 512421e10b8SBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 513421e10b8SBarry Smith } 514421e10b8SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 515421e10b8SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 516ccb205f8SBarry Smith CHKMEMQ; 517421e10b8SBarry Smith PetscFunctionReturn(0); 518421e10b8SBarry Smith } 519421e10b8SBarry Smith 520421e10b8SBarry Smith #undef __FUNCT__ 521290bbb0aSBarry Smith #define __FUNCT__ "MatMult_BlockMat_Symmetric" 522290bbb0aSBarry Smith PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y) 523290bbb0aSBarry Smith { 524290bbb0aSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 525290bbb0aSBarry Smith PetscErrorCode ierr; 526290bbb0aSBarry Smith PetscScalar *xx,*yy; 527d0f46423SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j; 528290bbb0aSBarry Smith Mat *aa; 529290bbb0aSBarry Smith 530290bbb0aSBarry Smith PetscFunctionBegin; 531290bbb0aSBarry Smith CHKMEMQ; 532290bbb0aSBarry Smith /* 533290bbb0aSBarry Smith Standard CSR multiply except each entry is a Mat 534290bbb0aSBarry Smith */ 535290bbb0aSBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 536290bbb0aSBarry Smith 537290bbb0aSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 538290bbb0aSBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 539290bbb0aSBarry Smith aj = bmat->j; 540290bbb0aSBarry Smith aa = bmat->a; 541290bbb0aSBarry Smith ii = bmat->i; 542290bbb0aSBarry Smith for (i=0; i<m; i++) { 543290bbb0aSBarry Smith jrow = ii[i]; 544290bbb0aSBarry Smith n = ii[i+1] - jrow; 545290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 546290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->middle,xx + bs*i);CHKERRQ(ierr); 547290bbb0aSBarry Smith /* if we ALWAYS required a diagonal entry then could remove this if test */ 548290bbb0aSBarry Smith if (aj[jrow] == i) { 549290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 550290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 551290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 552290bbb0aSBarry Smith jrow++; 553290bbb0aSBarry Smith n--; 554290bbb0aSBarry Smith } 555290bbb0aSBarry Smith for (j=0; j<n; j++) { 556290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); /* upper triangular part */ 557290bbb0aSBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 558290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 559290bbb0aSBarry Smith 560290bbb0aSBarry Smith ierr = VecPlaceArray(bmat->right,yy + bs*aj[jrow]);CHKERRQ(ierr); /* lower triangular part */ 561290bbb0aSBarry Smith ierr = MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);CHKERRQ(ierr); 562290bbb0aSBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 563290bbb0aSBarry Smith jrow++; 564290bbb0aSBarry Smith } 565290bbb0aSBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 566290bbb0aSBarry Smith ierr = VecResetArray(bmat->middle);CHKERRQ(ierr); 567290bbb0aSBarry Smith } 568290bbb0aSBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 569290bbb0aSBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 570290bbb0aSBarry Smith CHKMEMQ; 571290bbb0aSBarry Smith PetscFunctionReturn(0); 572290bbb0aSBarry Smith } 573290bbb0aSBarry Smith 574290bbb0aSBarry Smith #undef __FUNCT__ 575421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat" 576421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 577421e10b8SBarry Smith { 578421e10b8SBarry Smith PetscFunctionBegin; 579421e10b8SBarry Smith PetscFunctionReturn(0); 580421e10b8SBarry Smith } 581421e10b8SBarry Smith 582421e10b8SBarry Smith #undef __FUNCT__ 583421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat" 584421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 585421e10b8SBarry Smith { 586421e10b8SBarry Smith PetscFunctionBegin; 587421e10b8SBarry Smith PetscFunctionReturn(0); 588421e10b8SBarry Smith } 589421e10b8SBarry Smith 590421e10b8SBarry Smith #undef __FUNCT__ 591421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat" 592421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 593421e10b8SBarry Smith { 594421e10b8SBarry Smith PetscFunctionBegin; 595421e10b8SBarry Smith PetscFunctionReturn(0); 596421e10b8SBarry Smith } 597421e10b8SBarry Smith 598ccb205f8SBarry Smith /* 599ccb205f8SBarry Smith Adds diagonal pointers to sparse matrix structure. 600ccb205f8SBarry Smith */ 601ccb205f8SBarry Smith #undef __FUNCT__ 602ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat" 603ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) 604ccb205f8SBarry Smith { 605ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 606ccb205f8SBarry Smith PetscErrorCode ierr; 607d0f46423SBarry Smith PetscInt i,j,mbs = A->rmap->n/A->rmap->bs; 608ccb205f8SBarry Smith 609ccb205f8SBarry Smith PetscFunctionBegin; 610ccb205f8SBarry Smith if (!a->diag) { 611ccb205f8SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 612ccb205f8SBarry Smith } 613ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 614ccb205f8SBarry Smith a->diag[i] = a->i[i+1]; 615ccb205f8SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 616ccb205f8SBarry Smith if (a->j[j] == i) { 617ccb205f8SBarry Smith a->diag[i] = j; 618ccb205f8SBarry Smith break; 619ccb205f8SBarry Smith } 620ccb205f8SBarry Smith } 621ccb205f8SBarry Smith } 622ccb205f8SBarry Smith PetscFunctionReturn(0); 623ccb205f8SBarry Smith } 624ccb205f8SBarry Smith 625ccb205f8SBarry Smith #undef __FUNCT__ 626ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat" 6274aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 628ccb205f8SBarry Smith { 629ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 630ccb205f8SBarry Smith Mat_SeqAIJ *c; 631ccb205f8SBarry Smith PetscErrorCode ierr; 632ccb205f8SBarry Smith PetscInt i,k,first,step,lensi,nrows,ncols; 633d2a212eaSBarry Smith PetscInt *j_new,*i_new,*aj = a->j,*ailen = a->ilen; 6341620fd73SBarry Smith PetscScalar *a_new; 635ccb205f8SBarry Smith Mat C,*aa = a->a; 636ccb205f8SBarry Smith PetscTruth stride,equal; 637ccb205f8SBarry Smith 638ccb205f8SBarry Smith PetscFunctionBegin; 639ccb205f8SBarry Smith ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr); 640e32f2f54SBarry Smith if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for idential column and row indices"); 641ccb205f8SBarry Smith ierr = ISStride(iscol,&stride);CHKERRQ(ierr); 642e32f2f54SBarry Smith if (!stride) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices"); 643ccb205f8SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 644e32f2f54SBarry Smith if (step != A->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block"); 645ccb205f8SBarry Smith 646ccb205f8SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 647ccb205f8SBarry Smith ncols = nrows; 648ccb205f8SBarry Smith 649ccb205f8SBarry Smith /* create submatrix */ 650ccb205f8SBarry Smith if (scall == MAT_REUSE_MATRIX) { 651ccb205f8SBarry Smith PetscInt n_cols,n_rows; 652ccb205f8SBarry Smith C = *B; 653ccb205f8SBarry Smith ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr); 654e32f2f54SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 655ccb205f8SBarry Smith ierr = MatZeroEntries(C);CHKERRQ(ierr); 656ccb205f8SBarry Smith } else { 6577adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 658ccb205f8SBarry Smith ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 6596e63c7a1SBarry Smith if (A->symmetric) { 6606e63c7a1SBarry Smith ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 6616e63c7a1SBarry Smith } else { 662ccb205f8SBarry Smith ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr); 6636e63c7a1SBarry Smith } 6646e63c7a1SBarry Smith ierr = MatSeqAIJSetPreallocation(C,0,ailen);CHKERRQ(ierr); 6656e63c7a1SBarry Smith ierr = MatSeqSBAIJSetPreallocation(C,1,0,ailen);CHKERRQ(ierr); 666ccb205f8SBarry Smith } 667ccb205f8SBarry Smith c = (Mat_SeqAIJ*)C->data; 668ccb205f8SBarry Smith 669ccb205f8SBarry Smith /* loop over rows inserting into submatrix */ 670ccb205f8SBarry Smith a_new = c->a; 671ccb205f8SBarry Smith j_new = c->j; 672ccb205f8SBarry Smith i_new = c->i; 673ccb205f8SBarry Smith 674ccb205f8SBarry Smith for (i=0; i<nrows; i++) { 675ccb205f8SBarry Smith lensi = ailen[i]; 676ccb205f8SBarry Smith for (k=0; k<lensi; k++) { 677ccb205f8SBarry Smith *j_new++ = *aj++; 6781620fd73SBarry Smith ierr = MatGetValue(*aa++,first,first,a_new++);CHKERRQ(ierr); 679ccb205f8SBarry Smith } 680ccb205f8SBarry Smith i_new[i+1] = i_new[i] + lensi; 681ccb205f8SBarry Smith c->ilen[i] = lensi; 682ccb205f8SBarry Smith } 683ccb205f8SBarry Smith 684ccb205f8SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 685ccb205f8SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 686ccb205f8SBarry Smith *B = C; 687ccb205f8SBarry Smith PetscFunctionReturn(0); 688ccb205f8SBarry Smith } 689ccb205f8SBarry Smith 690421e10b8SBarry Smith #undef __FUNCT__ 691421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat" 692421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 693421e10b8SBarry Smith { 694421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 695421e10b8SBarry Smith PetscErrorCode ierr; 696421e10b8SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 697ccb205f8SBarry Smith PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0; 698421e10b8SBarry Smith Mat *aa = a->a,*ap; 699421e10b8SBarry Smith 700421e10b8SBarry Smith PetscFunctionBegin; 701421e10b8SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 702421e10b8SBarry Smith 703421e10b8SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 704421e10b8SBarry Smith for (i=1; i<m; i++) { 705421e10b8SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 706421e10b8SBarry Smith fshift += imax[i-1] - ailen[i-1]; 707421e10b8SBarry Smith rmax = PetscMax(rmax,ailen[i]); 708421e10b8SBarry Smith if (fshift) { 709421e10b8SBarry Smith ip = aj + ai[i] ; 710421e10b8SBarry Smith ap = aa + ai[i] ; 711421e10b8SBarry Smith N = ailen[i]; 712421e10b8SBarry Smith for (j=0; j<N; j++) { 713421e10b8SBarry Smith ip[j-fshift] = ip[j]; 714421e10b8SBarry Smith ap[j-fshift] = ap[j]; 715421e10b8SBarry Smith } 716421e10b8SBarry Smith } 717421e10b8SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 718421e10b8SBarry Smith } 719421e10b8SBarry Smith if (m) { 720421e10b8SBarry Smith fshift += imax[m-1] - ailen[m-1]; 721421e10b8SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 722421e10b8SBarry Smith } 723421e10b8SBarry Smith /* reset ilen and imax for each row */ 724421e10b8SBarry Smith for (i=0; i<m; i++) { 725421e10b8SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 726421e10b8SBarry Smith } 727421e10b8SBarry Smith a->nz = ai[m]; 728ccb205f8SBarry Smith for (i=0; i<a->nz; i++) { 729ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG) 730e32f2f54SBarry Smith if (!aa[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz); 731ccb205f8SBarry Smith #endif 732ccb205f8SBarry Smith ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 733ccb205f8SBarry Smith ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 734ccb205f8SBarry Smith } 735ccb205f8SBarry Smith CHKMEMQ; 736d0f46423SBarry 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); 737421e10b8SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); 738421e10b8SBarry Smith ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); 7398e58a170SBarry Smith A->info.mallocs += a->reallocs; 740421e10b8SBarry Smith a->reallocs = 0; 741421e10b8SBarry Smith A->info.nz_unneeded = (double)fshift; 742421e10b8SBarry Smith a->rmax = rmax; 743421e10b8SBarry Smith 744421e10b8SBarry Smith A->same_nonzero = PETSC_TRUE; 745ccb205f8SBarry Smith ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr); 746421e10b8SBarry Smith PetscFunctionReturn(0); 747421e10b8SBarry Smith } 748421e10b8SBarry Smith 749290bbb0aSBarry Smith #undef __FUNCT__ 750290bbb0aSBarry Smith #define __FUNCT__ "MatSetOption_BlockMat" 7514e0d8c25SBarry Smith PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscTruth flg) 752290bbb0aSBarry Smith { 753290bbb0aSBarry Smith PetscFunctionBegin; 7544e0d8c25SBarry Smith if (opt == MAT_SYMMETRIC && flg) { 75541f059aeSBarry Smith A->ops->sor = MatSOR_BlockMat_Symmetric; 756290bbb0aSBarry Smith A->ops->mult = MatMult_BlockMat_Symmetric; 757290bbb0aSBarry Smith } else { 758290bbb0aSBarry Smith PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]); 759290bbb0aSBarry Smith } 760290bbb0aSBarry Smith PetscFunctionReturn(0); 761290bbb0aSBarry Smith } 762290bbb0aSBarry Smith 763290bbb0aSBarry Smith 764421e10b8SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 765421e10b8SBarry Smith 0, 766421e10b8SBarry Smith 0, 767421e10b8SBarry Smith MatMult_BlockMat, 768421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat, 769421e10b8SBarry Smith MatMultTranspose_BlockMat, 770421e10b8SBarry Smith MatMultTransposeAdd_BlockMat, 771421e10b8SBarry Smith 0, 772421e10b8SBarry Smith 0, 773421e10b8SBarry Smith 0, 774421e10b8SBarry Smith /*10*/ 0, 775421e10b8SBarry Smith 0, 776421e10b8SBarry Smith 0, 77741f059aeSBarry Smith MatSOR_BlockMat, 778421e10b8SBarry Smith 0, 779421e10b8SBarry Smith /*15*/ 0, 780421e10b8SBarry Smith 0, 781421e10b8SBarry Smith 0, 782421e10b8SBarry Smith 0, 783421e10b8SBarry Smith 0, 784421e10b8SBarry Smith /*20*/ 0, 785421e10b8SBarry Smith MatAssemblyEnd_BlockMat, 786290bbb0aSBarry Smith MatSetOption_BlockMat, 787421e10b8SBarry Smith 0, 788d519adbfSMatthew Knepley /*24*/ 0, 789421e10b8SBarry Smith 0, 790421e10b8SBarry Smith 0, 791421e10b8SBarry Smith 0, 792421e10b8SBarry Smith 0, 793d519adbfSMatthew Knepley /*29*/ 0, 794421e10b8SBarry Smith 0, 795421e10b8SBarry Smith 0, 796421e10b8SBarry Smith 0, 797421e10b8SBarry Smith 0, 798d519adbfSMatthew Knepley /*34*/ 0, 799421e10b8SBarry Smith 0, 800421e10b8SBarry Smith 0, 801421e10b8SBarry Smith 0, 802421e10b8SBarry Smith 0, 803d519adbfSMatthew Knepley /*39*/ 0, 804421e10b8SBarry Smith 0, 805421e10b8SBarry Smith 0, 806421e10b8SBarry Smith 0, 807421e10b8SBarry Smith 0, 808d519adbfSMatthew Knepley /*44*/ 0, 809421e10b8SBarry Smith 0, 810421e10b8SBarry Smith 0, 811421e10b8SBarry Smith 0, 812421e10b8SBarry Smith 0, 813d519adbfSMatthew Knepley /*49*/ 0, 814421e10b8SBarry Smith 0, 815421e10b8SBarry Smith 0, 816421e10b8SBarry Smith 0, 817421e10b8SBarry Smith 0, 818d519adbfSMatthew Knepley /*54*/ 0, 819421e10b8SBarry Smith 0, 820421e10b8SBarry Smith 0, 821421e10b8SBarry Smith 0, 822421e10b8SBarry Smith 0, 823d519adbfSMatthew Knepley /*59*/ MatGetSubMatrix_BlockMat, 824421e10b8SBarry Smith MatDestroy_BlockMat, 825ccb205f8SBarry Smith MatView_BlockMat, 826421e10b8SBarry Smith 0, 827421e10b8SBarry Smith 0, 828d519adbfSMatthew Knepley /*64*/ 0, 829421e10b8SBarry Smith 0, 830421e10b8SBarry Smith 0, 831421e10b8SBarry Smith 0, 832421e10b8SBarry Smith 0, 833d519adbfSMatthew Knepley /*69*/ 0, 834421e10b8SBarry Smith 0, 835421e10b8SBarry Smith 0, 836421e10b8SBarry Smith 0, 837421e10b8SBarry Smith 0, 838d519adbfSMatthew Knepley /*74*/ 0, 839421e10b8SBarry Smith 0, 840421e10b8SBarry Smith 0, 841421e10b8SBarry Smith 0, 842421e10b8SBarry Smith 0, 843d519adbfSMatthew Knepley /*79*/ 0, 844421e10b8SBarry Smith 0, 845421e10b8SBarry Smith 0, 846421e10b8SBarry Smith 0, 847421e10b8SBarry Smith MatLoad_BlockMat, 848d519adbfSMatthew Knepley /*84*/ 0, 849421e10b8SBarry Smith 0, 850421e10b8SBarry Smith 0, 851421e10b8SBarry Smith 0, 852421e10b8SBarry Smith 0, 853d519adbfSMatthew Knepley /*89*/ 0, 854421e10b8SBarry Smith 0, 855421e10b8SBarry Smith 0, 856421e10b8SBarry Smith 0, 857421e10b8SBarry Smith 0, 858d519adbfSMatthew Knepley /*94*/ 0, 859421e10b8SBarry Smith 0, 860421e10b8SBarry Smith 0, 861421e10b8SBarry Smith 0}; 862421e10b8SBarry Smith 8638cdf2d9bSBarry Smith #undef __FUNCT__ 8648cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation" 8658cdf2d9bSBarry Smith /*@C 8668cdf2d9bSBarry Smith MatBlockMatSetPreallocation - For good matrix assembly performance 8678cdf2d9bSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 8688cdf2d9bSBarry Smith (or the array nnz). By setting these parameters accurately, performance 8698cdf2d9bSBarry Smith during matrix assembly can be increased by more than a factor of 50. 8708cdf2d9bSBarry Smith 8718cdf2d9bSBarry Smith Collective on MPI_Comm 8728cdf2d9bSBarry Smith 8738cdf2d9bSBarry Smith Input Parameters: 8748cdf2d9bSBarry Smith + B - The matrix 8758cdf2d9bSBarry Smith . bs - size of each block in matrix 8768cdf2d9bSBarry Smith . nz - number of nonzeros per block row (same for all rows) 8778cdf2d9bSBarry Smith - nnz - array containing the number of nonzeros in the various block rows 8788cdf2d9bSBarry Smith (possibly different for each row) or PETSC_NULL 8798cdf2d9bSBarry Smith 8808cdf2d9bSBarry Smith Notes: 8818cdf2d9bSBarry Smith If nnz is given then nz is ignored 8828cdf2d9bSBarry Smith 8838cdf2d9bSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 8848cdf2d9bSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 8858cdf2d9bSBarry Smith allocation. For large problems you MUST preallocate memory or you 8868cdf2d9bSBarry Smith will get TERRIBLE performance, see the users' manual chapter on matrices. 8878cdf2d9bSBarry Smith 8888cdf2d9bSBarry Smith Level: intermediate 8898cdf2d9bSBarry Smith 8908cdf2d9bSBarry Smith .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues() 8918cdf2d9bSBarry Smith 8928cdf2d9bSBarry Smith @*/ 8938cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 8948cdf2d9bSBarry Smith { 8958cdf2d9bSBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 8968cdf2d9bSBarry Smith 8978cdf2d9bSBarry Smith PetscFunctionBegin; 8988cdf2d9bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)B,"MatBlockMatSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 8998cdf2d9bSBarry Smith if (f) { 9008cdf2d9bSBarry Smith ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 9018cdf2d9bSBarry Smith } 9028cdf2d9bSBarry Smith PetscFunctionReturn(0); 9038cdf2d9bSBarry Smith } 9048cdf2d9bSBarry Smith 9058cdf2d9bSBarry Smith EXTERN_C_BEGIN 9068cdf2d9bSBarry Smith #undef __FUNCT__ 9078cdf2d9bSBarry Smith #define __FUNCT__ "MatBlockMatSetPreallocation_BlockMat" 9088cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz) 9098cdf2d9bSBarry Smith { 9108cdf2d9bSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 9118cdf2d9bSBarry Smith PetscErrorCode ierr; 9128cdf2d9bSBarry Smith PetscInt i; 9138cdf2d9bSBarry Smith 9148cdf2d9bSBarry Smith PetscFunctionBegin; 91565e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Block size given %D must be great than zero",bs); 916e32f2f54SBarry Smith if (A->rmap->n % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap->n); 917e32f2f54SBarry Smith if (A->cmap->n % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap->n); 9188cdf2d9bSBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 919e32f2f54SBarry Smith if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 9208cdf2d9bSBarry Smith if (nnz) { 921d0f46423SBarry Smith for (i=0; i<A->rmap->n/bs; i++) { 922e32f2f54SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 923e32f2f54SBarry Smith if (nnz[i] > A->cmap->n/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],A->cmap->n/bs); 9248cdf2d9bSBarry Smith } 9258cdf2d9bSBarry Smith } 926d0f46423SBarry Smith A->rmap->bs = A->cmap->bs = bs; 927d0f46423SBarry Smith bmat->mbs = A->rmap->n/bs; 9288cdf2d9bSBarry Smith 9298cdf2d9bSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr); 9308cdf2d9bSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->middle);CHKERRQ(ierr); 9318cdf2d9bSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr); 9328cdf2d9bSBarry Smith 9332ee49352SLisandro Dalcin if (!bmat->imax) { 934d0f46423SBarry Smith ierr = PetscMalloc2(A->rmap->n,PetscInt,&bmat->imax,A->rmap->n,PetscInt,&bmat->ilen);CHKERRQ(ierr); 935d0f46423SBarry Smith ierr = PetscLogObjectMemory(A,2*A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 9362ee49352SLisandro Dalcin } 9378cdf2d9bSBarry Smith if (nnz) { 9388cdf2d9bSBarry Smith nz = 0; 939d0f46423SBarry Smith for (i=0; i<A->rmap->n/A->rmap->bs; i++) { 9408cdf2d9bSBarry Smith bmat->imax[i] = nnz[i]; 9418cdf2d9bSBarry Smith nz += nnz[i]; 9428cdf2d9bSBarry Smith } 9438cdf2d9bSBarry Smith } else { 944e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation"); 9458cdf2d9bSBarry Smith } 9468cdf2d9bSBarry Smith 9478cdf2d9bSBarry Smith /* bmat->ilen will count nonzeros in each row so far. */ 9488cdf2d9bSBarry Smith for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;} 9498cdf2d9bSBarry Smith 9508cdf2d9bSBarry Smith /* allocate the matrix space */ 9512ee49352SLisandro Dalcin ierr = MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 952d0f46423SBarry Smith ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap->n+1,PetscInt,&bmat->i);CHKERRQ(ierr); 953d0f46423SBarry Smith ierr = PetscLogObjectMemory(A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 9548cdf2d9bSBarry Smith bmat->i[0] = 0; 9558cdf2d9bSBarry Smith for (i=1; i<bmat->mbs+1; i++) { 9568cdf2d9bSBarry Smith bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1]; 9578cdf2d9bSBarry Smith } 9588cdf2d9bSBarry Smith bmat->singlemalloc = PETSC_TRUE; 9598cdf2d9bSBarry Smith bmat->free_a = PETSC_TRUE; 9608cdf2d9bSBarry Smith bmat->free_ij = PETSC_TRUE; 9618cdf2d9bSBarry Smith 9628cdf2d9bSBarry Smith bmat->nz = 0; 9638cdf2d9bSBarry Smith bmat->maxnz = nz; 9648cdf2d9bSBarry Smith A->info.nz_unneeded = (double)bmat->maxnz; 9658cdf2d9bSBarry Smith 9668cdf2d9bSBarry Smith PetscFunctionReturn(0); 9678cdf2d9bSBarry Smith } 9688cdf2d9bSBarry Smith EXTERN_C_END 9698cdf2d9bSBarry Smith 970421e10b8SBarry Smith /*MC 971421e10b8SBarry Smith MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 972421e10b8SBarry Smith consisting of (usually) sparse blocks. 973421e10b8SBarry Smith 974421e10b8SBarry Smith Level: advanced 975421e10b8SBarry Smith 976421e10b8SBarry Smith .seealso: MatCreateBlockMat() 977421e10b8SBarry Smith 978421e10b8SBarry Smith M*/ 979421e10b8SBarry Smith 980421e10b8SBarry Smith EXTERN_C_BEGIN 981421e10b8SBarry Smith #undef __FUNCT__ 982421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat" 983421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A) 984421e10b8SBarry Smith { 985421e10b8SBarry Smith Mat_BlockMat *b; 986421e10b8SBarry Smith PetscErrorCode ierr; 987421e10b8SBarry Smith 988421e10b8SBarry Smith PetscFunctionBegin; 98938f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_BlockMat,&b);CHKERRQ(ierr); 990421e10b8SBarry Smith A->data = (void*)b; 99138f2d2fdSLisandro Dalcin ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 992421e10b8SBarry Smith 99326283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 99426283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 99526283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 99626283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 997421e10b8SBarry Smith 998421e10b8SBarry Smith A->assembled = PETSC_TRUE; 999421e10b8SBarry Smith A->preallocated = PETSC_FALSE; 1000421e10b8SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr); 1001290bbb0aSBarry Smith 10028cdf2d9bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatBlockMatSetPreallocation_C", 10038cdf2d9bSBarry Smith "MatBlockMatSetPreallocation_BlockMat", 10048cdf2d9bSBarry Smith MatBlockMatSetPreallocation_BlockMat);CHKERRQ(ierr); 10058cdf2d9bSBarry Smith 1006421e10b8SBarry Smith PetscFunctionReturn(0); 1007421e10b8SBarry Smith } 1008421e10b8SBarry Smith EXTERN_C_END 1009421e10b8SBarry Smith 1010421e10b8SBarry Smith #undef __FUNCT__ 1011421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat" 1012421e10b8SBarry Smith /*@C 1013421e10b8SBarry Smith MatCreateBlockMat - Creates a new matrix based sparse Mat storage 1014421e10b8SBarry Smith 1015421e10b8SBarry Smith Collective on MPI_Comm 1016421e10b8SBarry Smith 1017421e10b8SBarry Smith Input Parameters: 1018421e10b8SBarry Smith + comm - MPI communicator 1019421e10b8SBarry Smith . m - number of rows 1020421e10b8SBarry Smith . n - number of columns 10218cdf2d9bSBarry Smith . bs - size of each submatrix 10228cdf2d9bSBarry Smith . nz - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known) 10238cdf2d9bSBarry Smith - nnz - expected number of nonzers per block row if known (use PETSC_NULL otherwise) 1024421e10b8SBarry Smith 1025421e10b8SBarry Smith 1026421e10b8SBarry Smith Output Parameter: 1027421e10b8SBarry Smith . A - the matrix 1028421e10b8SBarry Smith 1029421e10b8SBarry Smith Level: intermediate 1030421e10b8SBarry Smith 1031421e10b8SBarry Smith PETSc requires that matrices and vectors being used for certain 1032421e10b8SBarry Smith operations are partitioned accordingly. For example, when 1033421e10b8SBarry Smith creating a bmat matrix, A, that supports parallel matrix-vector 1034421e10b8SBarry Smith products using MatMult(A,x,y) the user should set the number 1035421e10b8SBarry Smith of local matrix rows to be the number of local elements of the 1036421e10b8SBarry Smith corresponding result vector, y. Note that this is information is 1037421e10b8SBarry Smith required for use of the matrix interface routines, even though 1038421e10b8SBarry Smith the bmat matrix may not actually be physically partitioned. 1039421e10b8SBarry Smith For example, 1040421e10b8SBarry Smith 1041421e10b8SBarry Smith .keywords: matrix, bmat, create 1042421e10b8SBarry Smith 1043421e10b8SBarry Smith .seealso: MATBLOCKMAT 1044421e10b8SBarry Smith @*/ 10458cdf2d9bSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A) 1046421e10b8SBarry Smith { 1047421e10b8SBarry Smith PetscErrorCode ierr; 1048421e10b8SBarry Smith 1049421e10b8SBarry Smith PetscFunctionBegin; 1050421e10b8SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1051421e10b8SBarry Smith ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1052421e10b8SBarry Smith ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr); 10538cdf2d9bSBarry Smith ierr = MatBlockMatSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 1054421e10b8SBarry Smith PetscFunctionReturn(0); 1055421e10b8SBarry Smith } 1056421e10b8SBarry Smith 1057421e10b8SBarry Smith 1058421e10b8SBarry Smith 1059