1421e10b8SBarry Smith #define PETSCMAT_DLL 2421e10b8SBarry Smith 3421e10b8SBarry Smith /* 4421e10b8SBarry Smith This provides a matrix that consists of Mats 5421e10b8SBarry Smith */ 6421e10b8SBarry Smith 7421e10b8SBarry Smith #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 8421e10b8SBarry Smith #include "src/mat/impls/baij/seq/baij.h" /* use the common AIJ data-structure */ 9*ccb205f8SBarry Smith #include "petscksp.h" 10421e10b8SBarry Smith 11421e10b8SBarry Smith #define CHUNKSIZE 15 12421e10b8SBarry Smith 13421e10b8SBarry Smith typedef struct { 14421e10b8SBarry Smith SEQAIJHEADER(Mat); 15421e10b8SBarry Smith SEQBAIJHEADER; 16*ccb205f8SBarry Smith Mat *diags; 17421e10b8SBarry Smith 18421e10b8SBarry Smith Vec left,right; /* dummy vectors to perform local parts of product */ 19421e10b8SBarry Smith } Mat_BlockMat; 20421e10b8SBarry Smith 21421e10b8SBarry Smith #undef __FUNCT__ 22*ccb205f8SBarry Smith #define __FUNCT__ "MatRelax_BlockMat" 23*ccb205f8SBarry Smith PetscErrorCode MatRelax_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 24*ccb205f8SBarry Smith { 25*ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 26*ccb205f8SBarry Smith PetscScalar *x; 27*ccb205f8SBarry Smith const Mat *v = a->a; 28*ccb205f8SBarry Smith const PetscScalar *b; 29*ccb205f8SBarry Smith PetscErrorCode ierr; 30*ccb205f8SBarry Smith PetscInt n = A->cmap.n,i,mbs = n/A->rmap.bs,j,bs = A->rmap.bs; 31*ccb205f8SBarry Smith const PetscInt *idx; 32*ccb205f8SBarry Smith IS row,col; 33*ccb205f8SBarry Smith MatFactorInfo info; 34*ccb205f8SBarry Smith Vec left = a->left,right = a->right; 35*ccb205f8SBarry Smith Mat *diag; 36*ccb205f8SBarry Smith 37*ccb205f8SBarry Smith PetscFunctionBegin; 38*ccb205f8SBarry Smith its = its*lits; 39*ccb205f8SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 40*ccb205f8SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat"); 41*ccb205f8SBarry Smith if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support yet for omega not equal to 1.0"); 42*ccb205f8SBarry Smith if (fshift) SETERRQ(PETSC_ERR_SUP,"No support yet for fshift"); 43*ccb205f8SBarry Smith 44*ccb205f8SBarry Smith if (!a->diags) { 45*ccb205f8SBarry Smith ierr = PetscMalloc(mbs*sizeof(Mat),&a->diags);CHKERRQ(ierr); 46*ccb205f8SBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 47*ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 48*ccb205f8SBarry Smith ierr = MatGetOrdering(a->a[a->diag[i]], MATORDERING_ND,&row,&col);CHKERRQ(ierr); 49*ccb205f8SBarry Smith ierr = MatLUFactorSymbolic(a->a[a->diag[i]],row,col,&info,a->diags+i);CHKERRQ(ierr); 50*ccb205f8SBarry Smith ierr = MatLUFactorNumeric(a->a[a->diag[i]],&info,a->diags+i);CHKERRQ(ierr); 51*ccb205f8SBarry Smith ierr = ISDestroy(row);CHKERRQ(ierr); 52*ccb205f8SBarry Smith ierr = ISDestroy(col);CHKERRQ(ierr); 53*ccb205f8SBarry Smith } 54*ccb205f8SBarry Smith } 55*ccb205f8SBarry Smith diag = a->diags; 56*ccb205f8SBarry Smith 57*ccb205f8SBarry Smith ierr = VecSet(xx,0.0);CHKERRQ(ierr); 58*ccb205f8SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 59*ccb205f8SBarry Smith ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 60*ccb205f8SBarry Smith 61*ccb205f8SBarry Smith /* need to add code for when initial guess is zero, see MatRelax_SeqAIJ */ 62*ccb205f8SBarry Smith while (its--) { 63*ccb205f8SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 64*ccb205f8SBarry Smith 65*ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 66*ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 67*ccb205f8SBarry Smith idx = a->j + a->i[i]; 68*ccb205f8SBarry Smith v = a->a + a->i[i]; 69*ccb205f8SBarry Smith 70*ccb205f8SBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 71*ccb205f8SBarry Smith for (j=0; j<n; j++) { 72*ccb205f8SBarry Smith if (idx[j] != i) { 73*ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 74*ccb205f8SBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 75*ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 76*ccb205f8SBarry Smith } 77*ccb205f8SBarry Smith } 78*ccb205f8SBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 79*ccb205f8SBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 80*ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 81*ccb205f8SBarry Smith 82*ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 83*ccb205f8SBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 84*ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 85*ccb205f8SBarry Smith } 86*ccb205f8SBarry Smith } 87*ccb205f8SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 88*ccb205f8SBarry Smith 89*ccb205f8SBarry Smith for (i=mbs-1; i>=0; i--) { 90*ccb205f8SBarry Smith n = a->i[i+1] - a->i[i]; 91*ccb205f8SBarry Smith idx = a->j + a->i[i]; 92*ccb205f8SBarry Smith v = a->a + a->i[i]; 93*ccb205f8SBarry Smith 94*ccb205f8SBarry Smith ierr = VecSet(left,0.0);CHKERRQ(ierr); 95*ccb205f8SBarry Smith for (j=0; j<n; j++) { 96*ccb205f8SBarry Smith if (idx[j] != i) { 97*ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + idx[j]*bs);CHKERRQ(ierr); 98*ccb205f8SBarry Smith ierr = MatMultAdd(v[j],right,left,left);CHKERRQ(ierr); 99*ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 100*ccb205f8SBarry Smith } 101*ccb205f8SBarry Smith } 102*ccb205f8SBarry Smith ierr = VecPlaceArray(right,b + i*bs);CHKERRQ(ierr); 103*ccb205f8SBarry Smith ierr = VecAYPX(left,-1.0,right);CHKERRQ(ierr); 104*ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 105*ccb205f8SBarry Smith 106*ccb205f8SBarry Smith ierr = VecPlaceArray(right,x + i*bs);CHKERRQ(ierr); 107*ccb205f8SBarry Smith ierr = MatSolve(diag[i],left,right);CHKERRQ(ierr); 108*ccb205f8SBarry Smith ierr = VecResetArray(right);CHKERRQ(ierr); 109*ccb205f8SBarry Smith 110*ccb205f8SBarry Smith } 111*ccb205f8SBarry Smith } 112*ccb205f8SBarry Smith } 113*ccb205f8SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 114*ccb205f8SBarry Smith ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 115*ccb205f8SBarry Smith PetscFunctionReturn(0); 116*ccb205f8SBarry Smith } 117*ccb205f8SBarry Smith 118*ccb205f8SBarry Smith #undef __FUNCT__ 119421e10b8SBarry Smith #define __FUNCT__ "MatSetValues_BlockMat" 120421e10b8SBarry Smith PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 121421e10b8SBarry Smith { 122421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 123421e10b8SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 124421e10b8SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 125421e10b8SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol; 126421e10b8SBarry Smith PetscErrorCode ierr; 127421e10b8SBarry Smith PetscInt ridx,cidx; 128421e10b8SBarry Smith PetscTruth roworiented=a->roworiented; 129421e10b8SBarry Smith MatScalar value; 130421e10b8SBarry Smith Mat *ap,*aa = a->a; 131421e10b8SBarry Smith 132421e10b8SBarry Smith PetscFunctionBegin; 133421e10b8SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 134421e10b8SBarry Smith row = im[k]; 135421e10b8SBarry Smith brow = row/bs; 136421e10b8SBarry Smith if (row < 0) continue; 137421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 138421e10b8SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 139421e10b8SBarry Smith #endif 140421e10b8SBarry Smith rp = aj + ai[brow]; 141421e10b8SBarry Smith ap = aa + ai[brow]; 142421e10b8SBarry Smith rmax = imax[brow]; 143421e10b8SBarry Smith nrow = ailen[brow]; 144421e10b8SBarry Smith low = 0; 145421e10b8SBarry Smith high = nrow; 146421e10b8SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 147421e10b8SBarry Smith if (in[l] < 0) continue; 148421e10b8SBarry Smith #if defined(PETSC_USE_DEBUG) 149421e10b8SBarry 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); 150421e10b8SBarry Smith #endif 151421e10b8SBarry Smith col = in[l]; bcol = col/bs; 152421e10b8SBarry Smith ridx = row % bs; cidx = col % bs; 153421e10b8SBarry Smith if (roworiented) { 154421e10b8SBarry Smith value = v[l + k*n]; 155421e10b8SBarry Smith } else { 156421e10b8SBarry Smith value = v[k + l*m]; 157421e10b8SBarry Smith } 158421e10b8SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 159421e10b8SBarry Smith lastcol = col; 160421e10b8SBarry Smith while (high-low > 7) { 161421e10b8SBarry Smith t = (low+high)/2; 162421e10b8SBarry Smith if (rp[t] > bcol) high = t; 163421e10b8SBarry Smith else low = t; 164421e10b8SBarry Smith } 165421e10b8SBarry Smith for (i=low; i<high; i++) { 166421e10b8SBarry Smith if (rp[i] > bcol) break; 167421e10b8SBarry Smith if (rp[i] == bcol) { 168*ccb205f8SBarry Smith /* printf("row %d col %d found i %d\n",brow,bcol,i);*/ 169421e10b8SBarry Smith goto noinsert1; 170421e10b8SBarry Smith } 171421e10b8SBarry Smith } 172421e10b8SBarry Smith if (nonew == 1) goto noinsert1; 173421e10b8SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 174421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat); 175421e10b8SBarry Smith N = nrow++ - 1; high++; 176421e10b8SBarry Smith /* shift up all the later entries in this row */ 177*ccb205f8SBarry Smith /* printf("N %d i %d\n",N,i);*/ 178421e10b8SBarry Smith for (ii=N; ii>=i; ii--) { 179421e10b8SBarry Smith rp[ii+1] = rp[ii]; 180421e10b8SBarry Smith ap[ii+1] = ap[ii]; 181421e10b8SBarry Smith } 182421e10b8SBarry Smith if (N>=i) ap[i] = 0; 183421e10b8SBarry Smith rp[i] = bcol; 184421e10b8SBarry Smith a->nz++; 185421e10b8SBarry Smith noinsert1:; 186421e10b8SBarry Smith if (!*(ap+i)) { 187*ccb205f8SBarry Smith /* printf("create matrix at i %d rw %d col %d\n",i,brow,bcol);*/ 188421e10b8SBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr); 189421e10b8SBarry Smith } 190*ccb205f8SBarry Smith /* printf("numerical value at i %d row %d col %d cidx %d ridx %d value %g\n",i,brow,bcol,cidx,ridx,value);*/ 191421e10b8SBarry Smith ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr); 192421e10b8SBarry Smith low = i; 193421e10b8SBarry Smith } 194*ccb205f8SBarry Smith /* printf("nrow for row %d %d\n",nrow,brow);*/ 195421e10b8SBarry Smith ailen[brow] = nrow; 196421e10b8SBarry Smith } 197421e10b8SBarry Smith A->same_nonzero = PETSC_FALSE; 198421e10b8SBarry Smith PetscFunctionReturn(0); 199421e10b8SBarry Smith } 200421e10b8SBarry Smith 201421e10b8SBarry Smith #undef __FUNCT__ 202421e10b8SBarry Smith #define __FUNCT__ "MatLoad_BlockMat" 203421e10b8SBarry Smith PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, MatType type,Mat *A) 204421e10b8SBarry Smith { 205421e10b8SBarry Smith PetscErrorCode ierr; 206421e10b8SBarry Smith Mat tmpA; 207421e10b8SBarry Smith PetscInt i,m,n,bs = 1,ncols; 208421e10b8SBarry Smith const PetscInt *cols; 209421e10b8SBarry Smith const PetscScalar *values; 210421e10b8SBarry Smith 211421e10b8SBarry Smith PetscFunctionBegin; 212421e10b8SBarry Smith ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr); 213421e10b8SBarry Smith 214421e10b8SBarry Smith ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr); 215421e10b8SBarry Smith ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix","Mat");CHKERRQ(ierr); 216421e10b8SBarry Smith ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr); 217421e10b8SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 218421e10b8SBarry Smith ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,A);CHKERRQ(ierr); 219421e10b8SBarry Smith for (i=0; i<m; i++) { 220421e10b8SBarry Smith ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 221421e10b8SBarry Smith ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr); 222*ccb205f8SBarry Smith ierr = MatRestoreRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 223421e10b8SBarry Smith } 224421e10b8SBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 225421e10b8SBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 226421e10b8SBarry Smith PetscFunctionReturn(0); 227421e10b8SBarry Smith } 228421e10b8SBarry Smith 229*ccb205f8SBarry Smith #undef __FUNCT__ 230*ccb205f8SBarry Smith #define __FUNCT__ "MatView_BlockMat" 231*ccb205f8SBarry Smith PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer) 232*ccb205f8SBarry Smith { 233*ccb205f8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 234*ccb205f8SBarry Smith PetscErrorCode ierr; 235*ccb205f8SBarry Smith const char *name; 236*ccb205f8SBarry Smith PetscViewerFormat format; 237*ccb205f8SBarry Smith 238*ccb205f8SBarry Smith PetscFunctionBegin; 239*ccb205f8SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 240*ccb205f8SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 241*ccb205f8SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 242*ccb205f8SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);CHKERRQ(ierr); 243*ccb205f8SBarry Smith } 244*ccb205f8SBarry Smith PetscFunctionReturn(0); 245*ccb205f8SBarry Smith } 246421e10b8SBarry Smith 247421e10b8SBarry Smith #undef __FUNCT__ 248421e10b8SBarry Smith #define __FUNCT__ "MatDestroy_BlockMat" 249421e10b8SBarry Smith PetscErrorCode MatDestroy_BlockMat(Mat mat) 250421e10b8SBarry Smith { 251421e10b8SBarry Smith PetscErrorCode ierr; 252421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 253*ccb205f8SBarry Smith PetscInt i; 254421e10b8SBarry Smith 255421e10b8SBarry Smith PetscFunctionBegin; 256421e10b8SBarry Smith if (bmat->right) { 257421e10b8SBarry Smith ierr = VecDestroy(bmat->right);CHKERRQ(ierr); 258421e10b8SBarry Smith } 259421e10b8SBarry Smith if (bmat->left) { 260421e10b8SBarry Smith ierr = VecDestroy(bmat->left);CHKERRQ(ierr); 261421e10b8SBarry Smith } 262*ccb205f8SBarry Smith if (bmat->diags) { 263*ccb205f8SBarry Smith for (i=0; i<mat->rmap.n/mat->rmap.bs; i++) { 264*ccb205f8SBarry Smith if (bmat->diags[i]) {ierr = MatDestroy(bmat->diags[i]);CHKERRQ(ierr);} 265*ccb205f8SBarry Smith } 266*ccb205f8SBarry Smith } 267*ccb205f8SBarry Smith if (bmat->a) { 268*ccb205f8SBarry Smith for (i=0; i<bmat->nz; i++) { 269*ccb205f8SBarry Smith if (bmat->a[i]) {ierr = MatDestroy(bmat->a[i]);CHKERRQ(ierr);} 270*ccb205f8SBarry Smith } 271*ccb205f8SBarry Smith } 272*ccb205f8SBarry Smith ierr = MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);CHKERRQ(ierr); 273421e10b8SBarry Smith ierr = PetscFree(bmat);CHKERRQ(ierr); 274421e10b8SBarry Smith PetscFunctionReturn(0); 275421e10b8SBarry Smith } 276421e10b8SBarry Smith 277421e10b8SBarry Smith #undef __FUNCT__ 278421e10b8SBarry Smith #define __FUNCT__ "MatMult_BlockMat" 279421e10b8SBarry Smith PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 280421e10b8SBarry Smith { 281421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 282421e10b8SBarry Smith PetscErrorCode ierr; 283421e10b8SBarry Smith PetscScalar *xx,*yy; 284421e10b8SBarry Smith PetscInt *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j; 285421e10b8SBarry Smith Mat *aa; 286421e10b8SBarry Smith 287421e10b8SBarry Smith PetscFunctionBegin; 288*ccb205f8SBarry Smith CHKMEMQ; 289421e10b8SBarry Smith /* 290421e10b8SBarry Smith Standard CSR multiply except each entry is a Mat 291421e10b8SBarry Smith */ 292421e10b8SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 293*ccb205f8SBarry Smith 294421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 295421e10b8SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 296421e10b8SBarry Smith aj = bmat->j; 297421e10b8SBarry Smith aa = bmat->a; 298421e10b8SBarry Smith ii = bmat->i; 299421e10b8SBarry Smith for (i=0; i<m; i++) { 300421e10b8SBarry Smith jrow = ii[i]; 301*ccb205f8SBarry Smith ierr = VecPlaceArray(bmat->left,yy + bs*i);CHKERRQ(ierr); 302421e10b8SBarry Smith n = ii[i+1] - jrow; 303421e10b8SBarry Smith ierr = VecSet(bmat->left,0.0);CHKERRQ(ierr); 304421e10b8SBarry Smith for (j=0; j<n; j++) { 305421e10b8SBarry Smith ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 306*ccb205f8SBarry Smith ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);CHKERRQ(ierr); 307421e10b8SBarry Smith ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 308421e10b8SBarry Smith jrow++; 309421e10b8SBarry Smith } 310421e10b8SBarry Smith ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 311421e10b8SBarry Smith } 312421e10b8SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 313421e10b8SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 314*ccb205f8SBarry Smith CHKMEMQ; 315421e10b8SBarry Smith PetscFunctionReturn(0); 316421e10b8SBarry Smith } 317421e10b8SBarry Smith 318421e10b8SBarry Smith #undef __FUNCT__ 319421e10b8SBarry Smith #define __FUNCT__ "MatMultAdd_BlockMat" 320421e10b8SBarry Smith PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 321421e10b8SBarry Smith { 322421e10b8SBarry Smith PetscFunctionBegin; 323421e10b8SBarry Smith PetscFunctionReturn(0); 324421e10b8SBarry Smith } 325421e10b8SBarry Smith 326421e10b8SBarry Smith #undef __FUNCT__ 327421e10b8SBarry Smith #define __FUNCT__ "MatMultTranspose_BlockMat" 328421e10b8SBarry Smith PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 329421e10b8SBarry Smith { 330421e10b8SBarry Smith PetscFunctionBegin; 331421e10b8SBarry Smith PetscFunctionReturn(0); 332421e10b8SBarry Smith } 333421e10b8SBarry Smith 334421e10b8SBarry Smith #undef __FUNCT__ 335421e10b8SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_BlockMat" 336421e10b8SBarry Smith PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 337421e10b8SBarry Smith { 338421e10b8SBarry Smith PetscFunctionBegin; 339421e10b8SBarry Smith PetscFunctionReturn(0); 340421e10b8SBarry Smith } 341421e10b8SBarry Smith 342421e10b8SBarry Smith #undef __FUNCT__ 343421e10b8SBarry Smith #define __FUNCT__ "MatSetBlockSize_BlockMat" 344421e10b8SBarry Smith PetscErrorCode MatSetBlockSize_BlockMat(Mat A,PetscInt bs) 345421e10b8SBarry Smith { 346421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 347421e10b8SBarry Smith PetscErrorCode ierr; 348*ccb205f8SBarry Smith PetscInt nz = 10,i; 349421e10b8SBarry Smith 350421e10b8SBarry Smith PetscFunctionBegin; 351*ccb205f8SBarry Smith if (A->rmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of rows %D",bs,A->rmap.n); 352*ccb205f8SBarry Smith if (A->cmap.n % bs) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Blocksize %D does not divide number of columns %D",bs,A->cmap.n); 353421e10b8SBarry Smith A->rmap.bs = A->cmap.bs = bs; 354421e10b8SBarry Smith 355*ccb205f8SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr); 356*ccb205f8SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);CHKERRQ(ierr); 357*ccb205f8SBarry Smith 358*ccb205f8SBarry Smith 359421e10b8SBarry Smith ierr = PetscMalloc2(A->rmap.n,PetscInt,&bmat->imax,A->rmap.n,PetscInt,&bmat->ilen);CHKERRQ(ierr); 360421e10b8SBarry Smith for (i=0; i<A->rmap.n; i++) bmat->imax[i] = nz; 361421e10b8SBarry Smith nz = nz*A->rmap.n; 362421e10b8SBarry Smith 363*ccb205f8SBarry Smith bmat->mbs = A->rmap.n/A->rmap.bs; 364421e10b8SBarry Smith 365421e10b8SBarry Smith /* bmat->ilen will count nonzeros in each row so far. */ 366*ccb205f8SBarry Smith for (i=0; i<bmat->mbs; i++) { bmat->ilen[i] = 0;} 367421e10b8SBarry Smith 368421e10b8SBarry Smith /* allocate the matrix space */ 369421e10b8SBarry Smith ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap.n+1,PetscInt,&bmat->i);CHKERRQ(ierr); 370421e10b8SBarry Smith bmat->i[0] = 0; 371*ccb205f8SBarry Smith for (i=1; i<bmat->mbs+1; i++) { 372421e10b8SBarry Smith bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1]; 373421e10b8SBarry Smith } 374421e10b8SBarry Smith bmat->singlemalloc = PETSC_TRUE; 375421e10b8SBarry Smith bmat->free_a = PETSC_TRUE; 376421e10b8SBarry Smith bmat->free_ij = PETSC_TRUE; 377421e10b8SBarry Smith 378421e10b8SBarry Smith bmat->nz = 0; 379421e10b8SBarry Smith bmat->maxnz = nz; 380421e10b8SBarry Smith A->info.nz_unneeded = (double)bmat->maxnz; 381421e10b8SBarry Smith 382421e10b8SBarry Smith PetscFunctionReturn(0); 383421e10b8SBarry Smith } 384421e10b8SBarry Smith 385*ccb205f8SBarry Smith /* 386*ccb205f8SBarry Smith Adds diagonal pointers to sparse matrix structure. 387*ccb205f8SBarry Smith */ 388*ccb205f8SBarry Smith #undef __FUNCT__ 389*ccb205f8SBarry Smith #define __FUNCT__ "MatMarkDiagonal_BlockMat" 390*ccb205f8SBarry Smith PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) 391*ccb205f8SBarry Smith { 392*ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 393*ccb205f8SBarry Smith PetscErrorCode ierr; 394*ccb205f8SBarry Smith PetscInt i,j,mbs = A->rmap.n/A->rmap.bs; 395*ccb205f8SBarry Smith 396*ccb205f8SBarry Smith PetscFunctionBegin; 397*ccb205f8SBarry Smith if (!a->diag) { 398*ccb205f8SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 399*ccb205f8SBarry Smith } 400*ccb205f8SBarry Smith for (i=0; i<mbs; i++) { 401*ccb205f8SBarry Smith a->diag[i] = a->i[i+1]; 402*ccb205f8SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 403*ccb205f8SBarry Smith if (a->j[j] == i) { 404*ccb205f8SBarry Smith a->diag[i] = j; 405*ccb205f8SBarry Smith break; 406*ccb205f8SBarry Smith } 407*ccb205f8SBarry Smith } 408*ccb205f8SBarry Smith } 409*ccb205f8SBarry Smith PetscFunctionReturn(0); 410*ccb205f8SBarry Smith } 411*ccb205f8SBarry Smith 412*ccb205f8SBarry Smith #undef __FUNCT__ 413*ccb205f8SBarry Smith #define __FUNCT__ "MatGetSubMatrix_BlockMat" 414*ccb205f8SBarry Smith PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B) 415*ccb205f8SBarry Smith { 416*ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 417*ccb205f8SBarry Smith Mat_SeqAIJ *c; 418*ccb205f8SBarry Smith PetscErrorCode ierr; 419*ccb205f8SBarry Smith PetscInt i,k,first,step,lensi,nrows,ncols; 420*ccb205f8SBarry Smith PetscInt *j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen; 421*ccb205f8SBarry Smith PetscScalar *a_new,value; 422*ccb205f8SBarry Smith Mat C,*aa = a->a; 423*ccb205f8SBarry Smith PetscTruth stride,equal; 424*ccb205f8SBarry Smith 425*ccb205f8SBarry Smith PetscFunctionBegin; 426*ccb205f8SBarry Smith ierr = ISEqual(isrow,iscol,&equal);CHKERRQ(ierr); 427*ccb205f8SBarry Smith if (!equal) SETERRQ(PETSC_ERR_SUP,"Only for idential column and row indices"); 428*ccb205f8SBarry Smith ierr = ISStride(iscol,&stride);CHKERRQ(ierr); 429*ccb205f8SBarry Smith if (!stride) SETERRQ(PETSC_ERR_SUP,"Only for stride indices"); 430*ccb205f8SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); 431*ccb205f8SBarry Smith if (step != A->rmap.bs) SETERRQ(PETSC_ERR_SUP,"Can only select one entry from each block"); 432*ccb205f8SBarry Smith 433*ccb205f8SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 434*ccb205f8SBarry Smith ncols = nrows; 435*ccb205f8SBarry Smith 436*ccb205f8SBarry Smith /* create submatrix */ 437*ccb205f8SBarry Smith if (scall == MAT_REUSE_MATRIX) { 438*ccb205f8SBarry Smith PetscInt n_cols,n_rows; 439*ccb205f8SBarry Smith C = *B; 440*ccb205f8SBarry Smith ierr = MatGetSize(C,&n_rows,&n_cols);CHKERRQ(ierr); 441*ccb205f8SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 442*ccb205f8SBarry Smith ierr = MatZeroEntries(C);CHKERRQ(ierr); 443*ccb205f8SBarry Smith } else { 444*ccb205f8SBarry Smith ierr = MatCreate(A->comm,&C);CHKERRQ(ierr); 445*ccb205f8SBarry Smith ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 446*ccb205f8SBarry Smith ierr = MatSetType(C,MATSEQAIJ);CHKERRQ(ierr); 447*ccb205f8SBarry Smith ierr = MatSeqAIJSetPreallocation_SeqAIJ(C,0,ailen);CHKERRQ(ierr); 448*ccb205f8SBarry Smith } 449*ccb205f8SBarry Smith c = (Mat_SeqAIJ*)C->data; 450*ccb205f8SBarry Smith 451*ccb205f8SBarry Smith /* loop over rows inserting into submatrix */ 452*ccb205f8SBarry Smith a_new = c->a; 453*ccb205f8SBarry Smith j_new = c->j; 454*ccb205f8SBarry Smith i_new = c->i; 455*ccb205f8SBarry Smith 456*ccb205f8SBarry Smith for (i=0; i<nrows; i++) { 457*ccb205f8SBarry Smith ii = ai[i]; 458*ccb205f8SBarry Smith lensi = ailen[i]; 459*ccb205f8SBarry Smith for (k=0; k<lensi; k++) { 460*ccb205f8SBarry Smith *j_new++ = *aj++; 461*ccb205f8SBarry Smith ierr = MatGetValue(*aa++,first,first,value);CHKERRQ(ierr); 462*ccb205f8SBarry Smith *a_new++ = value; 463*ccb205f8SBarry Smith } 464*ccb205f8SBarry Smith i_new[i+1] = i_new[i] + lensi; 465*ccb205f8SBarry Smith c->ilen[i] = lensi; 466*ccb205f8SBarry Smith } 467*ccb205f8SBarry Smith 468*ccb205f8SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 469*ccb205f8SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 470*ccb205f8SBarry Smith *B = C; 471*ccb205f8SBarry Smith PetscFunctionReturn(0); 472*ccb205f8SBarry Smith } 473*ccb205f8SBarry Smith 474421e10b8SBarry Smith #undef __FUNCT__ 475421e10b8SBarry Smith #define __FUNCT__ "MatAssemblyEnd_BlockMat" 476421e10b8SBarry Smith PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 477421e10b8SBarry Smith { 478421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat*)A->data; 479421e10b8SBarry Smith PetscErrorCode ierr; 480421e10b8SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 481*ccb205f8SBarry Smith PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0; 482421e10b8SBarry Smith Mat *aa = a->a,*ap; 483421e10b8SBarry Smith 484421e10b8SBarry Smith PetscFunctionBegin; 485421e10b8SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 486421e10b8SBarry Smith 487421e10b8SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 488421e10b8SBarry Smith for (i=1; i<m; i++) { 489421e10b8SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 490421e10b8SBarry Smith fshift += imax[i-1] - ailen[i-1]; 491421e10b8SBarry Smith rmax = PetscMax(rmax,ailen[i]); 492421e10b8SBarry Smith if (fshift) { 493421e10b8SBarry Smith ip = aj + ai[i] ; 494421e10b8SBarry Smith ap = aa + ai[i] ; 495421e10b8SBarry Smith N = ailen[i]; 496421e10b8SBarry Smith for (j=0; j<N; j++) { 497421e10b8SBarry Smith ip[j-fshift] = ip[j]; 498421e10b8SBarry Smith ap[j-fshift] = ap[j]; 499421e10b8SBarry Smith } 500421e10b8SBarry Smith } 501421e10b8SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 502421e10b8SBarry Smith } 503421e10b8SBarry Smith if (m) { 504421e10b8SBarry Smith fshift += imax[m-1] - ailen[m-1]; 505421e10b8SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 506421e10b8SBarry Smith } 507421e10b8SBarry Smith /* reset ilen and imax for each row */ 508421e10b8SBarry Smith for (i=0; i<m; i++) { 509421e10b8SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 510421e10b8SBarry Smith } 511421e10b8SBarry Smith a->nz = ai[m]; 512*ccb205f8SBarry Smith for (i=0; i<a->nz; i++) { 513*ccb205f8SBarry Smith #if defined(PETSC_USE_DEBUG) 514*ccb205f8SBarry Smith if (!aa[i]) SETERRQ3(PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz); 515*ccb205f8SBarry Smith #endif 516*ccb205f8SBarry Smith ierr = MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 517*ccb205f8SBarry Smith ierr = MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 518*ccb205f8SBarry Smith } 519*ccb205f8SBarry Smith CHKMEMQ; 520421e10b8SBarry Smith ierr = PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);CHKERRQ(ierr); 521421e10b8SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); 522421e10b8SBarry Smith ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); 523421e10b8SBarry Smith a->reallocs = 0; 524421e10b8SBarry Smith A->info.nz_unneeded = (double)fshift; 525421e10b8SBarry Smith a->rmax = rmax; 526421e10b8SBarry Smith 527421e10b8SBarry Smith A->same_nonzero = PETSC_TRUE; 528*ccb205f8SBarry Smith ierr = MatMarkDiagonal_BlockMat(A);CHKERRQ(ierr); 529421e10b8SBarry Smith PetscFunctionReturn(0); 530421e10b8SBarry Smith } 531421e10b8SBarry Smith 532421e10b8SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 533421e10b8SBarry Smith 0, 534421e10b8SBarry Smith 0, 535421e10b8SBarry Smith MatMult_BlockMat, 536421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat, 537421e10b8SBarry Smith MatMultTranspose_BlockMat, 538421e10b8SBarry Smith MatMultTransposeAdd_BlockMat, 539421e10b8SBarry Smith 0, 540421e10b8SBarry Smith 0, 541421e10b8SBarry Smith 0, 542421e10b8SBarry Smith /*10*/ 0, 543421e10b8SBarry Smith 0, 544421e10b8SBarry Smith 0, 545*ccb205f8SBarry Smith MatRelax_BlockMat, 546421e10b8SBarry Smith 0, 547421e10b8SBarry Smith /*15*/ 0, 548421e10b8SBarry Smith 0, 549421e10b8SBarry Smith 0, 550421e10b8SBarry Smith 0, 551421e10b8SBarry Smith 0, 552421e10b8SBarry Smith /*20*/ 0, 553421e10b8SBarry Smith MatAssemblyEnd_BlockMat, 554421e10b8SBarry Smith 0, 555421e10b8SBarry Smith 0, 556421e10b8SBarry Smith 0, 557421e10b8SBarry Smith /*25*/ 0, 558421e10b8SBarry Smith 0, 559421e10b8SBarry Smith 0, 560421e10b8SBarry Smith 0, 561421e10b8SBarry Smith 0, 562421e10b8SBarry Smith /*30*/ 0, 563421e10b8SBarry Smith 0, 564421e10b8SBarry Smith 0, 565421e10b8SBarry Smith 0, 566421e10b8SBarry Smith 0, 567421e10b8SBarry Smith /*35*/ 0, 568421e10b8SBarry Smith 0, 569421e10b8SBarry Smith 0, 570421e10b8SBarry Smith 0, 571421e10b8SBarry Smith 0, 572421e10b8SBarry Smith /*40*/ 0, 573421e10b8SBarry Smith 0, 574421e10b8SBarry Smith 0, 575421e10b8SBarry Smith 0, 576421e10b8SBarry Smith 0, 577421e10b8SBarry Smith /*45*/ 0, 578421e10b8SBarry Smith 0, 579421e10b8SBarry Smith 0, 580421e10b8SBarry Smith 0, 581421e10b8SBarry Smith 0, 582421e10b8SBarry Smith /*50*/ MatSetBlockSize_BlockMat, 583421e10b8SBarry Smith 0, 584421e10b8SBarry Smith 0, 585421e10b8SBarry Smith 0, 586421e10b8SBarry Smith 0, 587421e10b8SBarry Smith /*55*/ 0, 588421e10b8SBarry Smith 0, 589421e10b8SBarry Smith 0, 590421e10b8SBarry Smith 0, 591421e10b8SBarry Smith 0, 592*ccb205f8SBarry Smith /*60*/ MatGetSubMatrix_BlockMat, 593421e10b8SBarry Smith MatDestroy_BlockMat, 594*ccb205f8SBarry Smith MatView_BlockMat, 595421e10b8SBarry Smith 0, 596421e10b8SBarry Smith 0, 597421e10b8SBarry Smith /*65*/ 0, 598421e10b8SBarry Smith 0, 599421e10b8SBarry Smith 0, 600421e10b8SBarry Smith 0, 601421e10b8SBarry Smith 0, 602421e10b8SBarry Smith /*70*/ 0, 603421e10b8SBarry Smith 0, 604421e10b8SBarry Smith 0, 605421e10b8SBarry Smith 0, 606421e10b8SBarry Smith 0, 607421e10b8SBarry Smith /*75*/ 0, 608421e10b8SBarry Smith 0, 609421e10b8SBarry Smith 0, 610421e10b8SBarry Smith 0, 611421e10b8SBarry Smith 0, 612421e10b8SBarry Smith /*80*/ 0, 613421e10b8SBarry Smith 0, 614421e10b8SBarry Smith 0, 615421e10b8SBarry Smith 0, 616421e10b8SBarry Smith MatLoad_BlockMat, 617421e10b8SBarry Smith /*85*/ 0, 618421e10b8SBarry Smith 0, 619421e10b8SBarry Smith 0, 620421e10b8SBarry Smith 0, 621421e10b8SBarry Smith 0, 622421e10b8SBarry Smith /*90*/ 0, 623421e10b8SBarry Smith 0, 624421e10b8SBarry Smith 0, 625421e10b8SBarry Smith 0, 626421e10b8SBarry Smith 0, 627421e10b8SBarry Smith /*95*/ 0, 628421e10b8SBarry Smith 0, 629421e10b8SBarry Smith 0, 630421e10b8SBarry Smith 0}; 631421e10b8SBarry Smith 632421e10b8SBarry Smith /*MC 633421e10b8SBarry Smith MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 634421e10b8SBarry Smith consisting of (usually) sparse blocks. 635421e10b8SBarry Smith 636421e10b8SBarry Smith Level: advanced 637421e10b8SBarry Smith 638421e10b8SBarry Smith .seealso: MatCreateBlockMat() 639421e10b8SBarry Smith 640421e10b8SBarry Smith M*/ 641421e10b8SBarry Smith 642421e10b8SBarry Smith EXTERN_C_BEGIN 643421e10b8SBarry Smith #undef __FUNCT__ 644421e10b8SBarry Smith #define __FUNCT__ "MatCreate_BlockMat" 645421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A) 646421e10b8SBarry Smith { 647421e10b8SBarry Smith Mat_BlockMat *b; 648421e10b8SBarry Smith PetscErrorCode ierr; 649421e10b8SBarry Smith 650421e10b8SBarry Smith PetscFunctionBegin; 651421e10b8SBarry Smith ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 652421e10b8SBarry Smith ierr = PetscNew(Mat_BlockMat,&b);CHKERRQ(ierr); 653421e10b8SBarry Smith 654421e10b8SBarry Smith A->data = (void*)b; 655421e10b8SBarry Smith 656421e10b8SBarry Smith ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr); 657421e10b8SBarry Smith ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr); 658421e10b8SBarry Smith 659421e10b8SBarry Smith A->assembled = PETSC_TRUE; 660421e10b8SBarry Smith A->preallocated = PETSC_FALSE; 661421e10b8SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr); 662421e10b8SBarry Smith PetscFunctionReturn(0); 663421e10b8SBarry Smith } 664421e10b8SBarry Smith EXTERN_C_END 665421e10b8SBarry Smith 666421e10b8SBarry Smith #undef __FUNCT__ 667421e10b8SBarry Smith #define __FUNCT__ "MatCreateBlockMat" 668421e10b8SBarry Smith /*@C 669421e10b8SBarry Smith MatCreateBlockMat - Creates a new matrix based sparse Mat storage 670421e10b8SBarry Smith 671421e10b8SBarry Smith Collective on MPI_Comm 672421e10b8SBarry Smith 673421e10b8SBarry Smith Input Parameters: 674421e10b8SBarry Smith + comm - MPI communicator 675421e10b8SBarry Smith . m - number of rows 676421e10b8SBarry Smith . n - number of columns 677421e10b8SBarry Smith - bs - size of each submatrix 678421e10b8SBarry Smith 679421e10b8SBarry Smith 680421e10b8SBarry Smith Output Parameter: 681421e10b8SBarry Smith . A - the matrix 682421e10b8SBarry Smith 683421e10b8SBarry Smith Level: intermediate 684421e10b8SBarry Smith 685421e10b8SBarry Smith PETSc requires that matrices and vectors being used for certain 686421e10b8SBarry Smith operations are partitioned accordingly. For example, when 687421e10b8SBarry Smith creating a bmat matrix, A, that supports parallel matrix-vector 688421e10b8SBarry Smith products using MatMult(A,x,y) the user should set the number 689421e10b8SBarry Smith of local matrix rows to be the number of local elements of the 690421e10b8SBarry Smith corresponding result vector, y. Note that this is information is 691421e10b8SBarry Smith required for use of the matrix interface routines, even though 692421e10b8SBarry Smith the bmat matrix may not actually be physically partitioned. 693421e10b8SBarry Smith For example, 694421e10b8SBarry Smith 695421e10b8SBarry Smith .keywords: matrix, bmat, create 696421e10b8SBarry Smith 697421e10b8SBarry Smith .seealso: MATBLOCKMAT 698421e10b8SBarry Smith @*/ 699421e10b8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,Mat *A) 700421e10b8SBarry Smith { 701421e10b8SBarry Smith PetscErrorCode ierr; 702421e10b8SBarry Smith 703421e10b8SBarry Smith PetscFunctionBegin; 704421e10b8SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 705421e10b8SBarry Smith ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 706421e10b8SBarry Smith ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr); 707421e10b8SBarry Smith ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 708421e10b8SBarry Smith PetscFunctionReturn(0); 709421e10b8SBarry Smith } 710421e10b8SBarry Smith 711421e10b8SBarry Smith 712421e10b8SBarry Smith 713