147c6ae99SBarry Smith 2c6db04a5SJed Brown #include <private/daimpl.h> /*I "petscdmda.h" I*/ 3c6db04a5SJed Brown #include <petscmat.h> /*I "petscmat.h" I*/ 4c6db04a5SJed Brown #include <private/matimpl.h> 547c6ae99SBarry Smith 6e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring *); 7e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring *); 8e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring *); 9e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM,ISColoringType,ISColoring *); 1047c6ae99SBarry Smith 1147c6ae99SBarry Smith /* 1247c6ae99SBarry Smith For ghost i that may be negative or greater than the upper bound this 1347c6ae99SBarry Smith maps it into the 0:m-1 range using periodicity 1447c6ae99SBarry Smith */ 1547c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i:((i >= m) ? i-m:i)) 1647c6ae99SBarry Smith 1747c6ae99SBarry Smith #undef __FUNCT__ 18aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills_Private" 19aa219208SBarry Smith static PetscErrorCode DMDASetBlockFills_Private(PetscInt *dfill,PetscInt w,PetscInt **rfill) 2047c6ae99SBarry Smith { 2147c6ae99SBarry Smith PetscErrorCode ierr; 2247c6ae99SBarry Smith PetscInt i,j,nz,*fill; 2347c6ae99SBarry Smith 2447c6ae99SBarry Smith PetscFunctionBegin; 2547c6ae99SBarry Smith if (!dfill) PetscFunctionReturn(0); 2647c6ae99SBarry Smith 2747c6ae99SBarry Smith /* count number nonzeros */ 2847c6ae99SBarry Smith nz = 0; 2947c6ae99SBarry Smith for (i=0; i<w; i++) { 3047c6ae99SBarry Smith for (j=0; j<w; j++) { 3147c6ae99SBarry Smith if (dfill[w*i+j]) nz++; 3247c6ae99SBarry Smith } 3347c6ae99SBarry Smith } 3447c6ae99SBarry Smith ierr = PetscMalloc((nz + w + 1)*sizeof(PetscInt),&fill);CHKERRQ(ierr); 3547c6ae99SBarry Smith /* construct modified CSR storage of nonzero structure */ 3647c6ae99SBarry Smith nz = w + 1; 3747c6ae99SBarry Smith for (i=0; i<w; i++) { 3847c6ae99SBarry Smith fill[i] = nz; 3947c6ae99SBarry Smith for (j=0; j<w; j++) { 4047c6ae99SBarry Smith if (dfill[w*i+j]) { 4147c6ae99SBarry Smith fill[nz] = j; 4247c6ae99SBarry Smith nz++; 4347c6ae99SBarry Smith } 4447c6ae99SBarry Smith } 4547c6ae99SBarry Smith } 4647c6ae99SBarry Smith fill[w] = nz; 4747c6ae99SBarry Smith 4847c6ae99SBarry Smith *rfill = fill; 4947c6ae99SBarry Smith PetscFunctionReturn(0); 5047c6ae99SBarry Smith } 5147c6ae99SBarry Smith 5247c6ae99SBarry Smith #undef __FUNCT__ 53aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills" 5447c6ae99SBarry Smith /*@ 55aa219208SBarry Smith DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem 56950540a4SJed Brown of the matrix returned by DMCreateMatrix(). 5747c6ae99SBarry Smith 58aa219208SBarry Smith Logically Collective on DMDA 5947c6ae99SBarry Smith 6047c6ae99SBarry Smith Input Parameter: 6147c6ae99SBarry Smith + da - the distributed array 6247c6ae99SBarry Smith . dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block) 6347c6ae99SBarry Smith - ofill - the fill pattern in the off-diagonal blocks 6447c6ae99SBarry Smith 6547c6ae99SBarry Smith 6647c6ae99SBarry Smith Level: developer 6747c6ae99SBarry Smith 6847c6ae99SBarry Smith Notes: This only makes sense when you are doing multicomponent problems but using the 6947c6ae99SBarry Smith MPIAIJ matrix format 7047c6ae99SBarry Smith 7147c6ae99SBarry Smith The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries 7247c6ae99SBarry Smith representing coupling and 0 entries for missing coupling. For example 7347c6ae99SBarry Smith $ dfill[9] = {1, 0, 0, 7447c6ae99SBarry Smith $ 1, 1, 0, 7547c6ae99SBarry Smith $ 0, 1, 1} 7647c6ae99SBarry Smith means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with 7747c6ae99SBarry Smith itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the 7847c6ae99SBarry Smith diagonal block). 7947c6ae99SBarry Smith 80aa219208SBarry Smith DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then 8147c6ae99SBarry Smith can be represented in the dfill, ofill format 8247c6ae99SBarry Smith 8347c6ae99SBarry Smith Contributed by Glenn Hammond 8447c6ae99SBarry Smith 85*8ddb5d8bSBarry Smith .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly() 8647c6ae99SBarry Smith 8747c6ae99SBarry Smith @*/ 887087cfbeSBarry Smith PetscErrorCode DMDASetBlockFills(DM da,PetscInt *dfill,PetscInt *ofill) 8947c6ae99SBarry Smith { 9047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 9147c6ae99SBarry Smith PetscErrorCode ierr; 9247c6ae99SBarry Smith 9347c6ae99SBarry Smith PetscFunctionBegin; 94aa219208SBarry Smith ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr); 95aa219208SBarry Smith ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr); 9647c6ae99SBarry Smith PetscFunctionReturn(0); 9747c6ae99SBarry Smith } 9847c6ae99SBarry Smith 9947c6ae99SBarry Smith 10047c6ae99SBarry Smith #undef __FUNCT__ 101e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA" 102e727c939SJed Brown PetscErrorCode DMCreateColoring_DA(DM da,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 10347c6ae99SBarry Smith { 10447c6ae99SBarry Smith PetscErrorCode ierr; 10547c6ae99SBarry Smith PetscInt dim,m,n,p,nc; 1061321219cSEthan Coon DMDABoundaryType bx,by,bz; 10747c6ae99SBarry Smith MPI_Comm comm; 10847c6ae99SBarry Smith PetscMPIInt size; 10947c6ae99SBarry Smith PetscBool isBAIJ; 11047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 11147c6ae99SBarry Smith 11247c6ae99SBarry Smith PetscFunctionBegin; 11347c6ae99SBarry Smith /* 11447c6ae99SBarry Smith m 11547c6ae99SBarry Smith ------------------------------------------------------ 11647c6ae99SBarry Smith | | 11747c6ae99SBarry Smith | | 11847c6ae99SBarry Smith | ---------------------- | 11947c6ae99SBarry Smith | | | | 12047c6ae99SBarry Smith n | yn | | | 12147c6ae99SBarry Smith | | | | 12247c6ae99SBarry Smith | .--------------------- | 12347c6ae99SBarry Smith | (xs,ys) xn | 12447c6ae99SBarry Smith | . | 12547c6ae99SBarry Smith | (gxs,gys) | 12647c6ae99SBarry Smith | | 12747c6ae99SBarry Smith ----------------------------------------------------- 12847c6ae99SBarry Smith */ 12947c6ae99SBarry Smith 13047c6ae99SBarry Smith /* 13147c6ae99SBarry Smith nc - number of components per grid point 13247c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 13347c6ae99SBarry Smith 13447c6ae99SBarry Smith */ 1351321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);CHKERRQ(ierr); 13647c6ae99SBarry Smith 13747c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 13847c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 13947c6ae99SBarry Smith if (ctype == IS_COLORING_GHOSTED){ 14047c6ae99SBarry Smith if (size == 1) { 14147c6ae99SBarry Smith ctype = IS_COLORING_GLOBAL; 14247c6ae99SBarry Smith } else if (dim > 1){ 1431321219cSEthan Coon if ((m==1 && bx == DMDA_BOUNDARY_PERIODIC) || (n==1 && by == DMDA_BOUNDARY_PERIODIC) || (p==1 && bz == DMDA_BOUNDARY_PERIODIC)){ 14447c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"IS_COLORING_GHOSTED cannot be used for periodic boundary condition having both ends of the domain on the same process"); 14547c6ae99SBarry Smith } 14647c6ae99SBarry Smith } 14747c6ae99SBarry Smith } 14847c6ae99SBarry Smith 149aa219208SBarry Smith /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ 15047c6ae99SBarry Smith matrices is for the blocks, not the individual matrix elements */ 1514833614aSSatish Balay ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr); 1524833614aSSatish Balay if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);} 15347c6ae99SBarry Smith if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);} 15447c6ae99SBarry Smith if (isBAIJ) { 15547c6ae99SBarry Smith dd->w = 1; 15647c6ae99SBarry Smith dd->xs = dd->xs/nc; 15747c6ae99SBarry Smith dd->xe = dd->xe/nc; 15847c6ae99SBarry Smith dd->Xs = dd->Xs/nc; 15947c6ae99SBarry Smith dd->Xe = dd->Xe/nc; 16047c6ae99SBarry Smith } 16147c6ae99SBarry Smith 16247c6ae99SBarry Smith /* 163aa219208SBarry Smith We do not provide a getcoloring function in the DMDA operations because 164aa219208SBarry Smith the basic DMDA does not know about matrices. We think of DMDA as being more 16547c6ae99SBarry Smith more low-level then matrices. 16647c6ae99SBarry Smith */ 16747c6ae99SBarry Smith if (dim == 1) { 168e727c939SJed Brown ierr = DMCreateColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 16947c6ae99SBarry Smith } else if (dim == 2) { 170e727c939SJed Brown ierr = DMCreateColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 17147c6ae99SBarry Smith } else if (dim == 3) { 172e727c939SJed Brown ierr = DMCreateColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 17371cd77b2SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim); 17447c6ae99SBarry Smith if (isBAIJ) { 17547c6ae99SBarry Smith dd->w = nc; 17647c6ae99SBarry Smith dd->xs = dd->xs*nc; 17747c6ae99SBarry Smith dd->xe = dd->xe*nc; 17847c6ae99SBarry Smith dd->Xs = dd->Xs*nc; 17947c6ae99SBarry Smith dd->Xe = dd->Xe*nc; 18047c6ae99SBarry Smith } 18147c6ae99SBarry Smith PetscFunctionReturn(0); 18247c6ae99SBarry Smith } 18347c6ae99SBarry Smith 18447c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 18547c6ae99SBarry Smith 18647c6ae99SBarry Smith #undef __FUNCT__ 187e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_MPIAIJ" 188e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 18947c6ae99SBarry Smith { 19047c6ae99SBarry Smith PetscErrorCode ierr; 19147c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col; 19247c6ae99SBarry Smith PetscInt ncolors; 19347c6ae99SBarry Smith MPI_Comm comm; 1941321219cSEthan Coon DMDABoundaryType bx,by; 195aa219208SBarry Smith DMDAStencilType st; 19647c6ae99SBarry Smith ISColoringValue *colors; 19747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 19847c6ae99SBarry Smith 19947c6ae99SBarry Smith PetscFunctionBegin; 20047c6ae99SBarry Smith /* 20147c6ae99SBarry Smith nc - number of components per grid point 20247c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 20347c6ae99SBarry Smith 20447c6ae99SBarry Smith */ 2051321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 20647c6ae99SBarry Smith col = 2*s + 1; 207aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 208aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 20947c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 21047c6ae99SBarry Smith 21147c6ae99SBarry Smith /* special case as taught to us by Paul Hovland */ 212aa219208SBarry Smith if (st == DMDA_STENCIL_STAR && s == 1) { 213e727c939SJed Brown ierr = DMCreateColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 21447c6ae99SBarry Smith } else { 21547c6ae99SBarry Smith 2161321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){ 21747c6ae99SBarry Smith SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\ 21847c6ae99SBarry Smith by 2*stencil_width + 1 (%d)\n", m, col); 21947c6ae99SBarry Smith } 2201321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){ 22147c6ae99SBarry Smith SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\ 22247c6ae99SBarry Smith by 2*stencil_width + 1 (%d)\n", n, col); 22347c6ae99SBarry Smith } 22447c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 22547c6ae99SBarry Smith if (!dd->localcoloring) { 22647c6ae99SBarry Smith ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 22747c6ae99SBarry Smith ii = 0; 22847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 22947c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 23047c6ae99SBarry Smith for (k=0; k<nc; k++) { 23147c6ae99SBarry Smith colors[ii++] = k + nc*((i % col) + col*(j % col)); 23247c6ae99SBarry Smith } 23347c6ae99SBarry Smith } 23447c6ae99SBarry Smith } 23547c6ae99SBarry Smith ncolors = nc + nc*(col-1 + col*(col-1)); 23647c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr); 23747c6ae99SBarry Smith } 23847c6ae99SBarry Smith *coloring = dd->localcoloring; 23947c6ae99SBarry Smith } else if (ctype == IS_COLORING_GHOSTED) { 24047c6ae99SBarry Smith if (!dd->ghostedcoloring) { 24147c6ae99SBarry Smith ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 24247c6ae99SBarry Smith ii = 0; 24347c6ae99SBarry Smith for (j=gys; j<gys+gny; j++) { 24447c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 24547c6ae99SBarry Smith for (k=0; k<nc; k++) { 24647c6ae99SBarry Smith /* the complicated stuff is to handle periodic boundaries */ 24747c6ae99SBarry Smith colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col)); 24847c6ae99SBarry Smith } 24947c6ae99SBarry Smith } 25047c6ae99SBarry Smith } 25147c6ae99SBarry Smith ncolors = nc + nc*(col - 1 + col*(col-1)); 25247c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr); 25347c6ae99SBarry Smith /* PetscIntView(ncolors,(PetscInt *)colors,0); */ 25447c6ae99SBarry Smith 25547c6ae99SBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr); 25647c6ae99SBarry Smith } 25747c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 25847c6ae99SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 25947c6ae99SBarry Smith } 26047c6ae99SBarry Smith ierr = ISColoringReference(*coloring);CHKERRQ(ierr); 26147c6ae99SBarry Smith PetscFunctionReturn(0); 26247c6ae99SBarry Smith } 26347c6ae99SBarry Smith 26447c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 26547c6ae99SBarry Smith 26647c6ae99SBarry Smith #undef __FUNCT__ 267e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_3d_MPIAIJ" 268e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 26947c6ae99SBarry Smith { 27047c6ae99SBarry Smith PetscErrorCode ierr; 27147c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,gxs,gys,gnx,gny,m,n,p,dim,s,k,nc,col,zs,gzs,ii,l,nz,gnz,M,N,P; 27247c6ae99SBarry Smith PetscInt ncolors; 27347c6ae99SBarry Smith MPI_Comm comm; 2741321219cSEthan Coon DMDABoundaryType bx,by,bz; 275aa219208SBarry Smith DMDAStencilType st; 27647c6ae99SBarry Smith ISColoringValue *colors; 27747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 27847c6ae99SBarry Smith 27947c6ae99SBarry Smith PetscFunctionBegin; 28047c6ae99SBarry Smith /* 28147c6ae99SBarry Smith nc - number of components per grid point 28247c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 28347c6ae99SBarry Smith 28447c6ae99SBarry Smith */ 2851321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 28647c6ae99SBarry Smith col = 2*s + 1; 2871321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){ 28847c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\ 28947c6ae99SBarry Smith by 2*stencil_width + 1\n"); 29047c6ae99SBarry Smith } 2911321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){ 29247c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\ 29347c6ae99SBarry Smith by 2*stencil_width + 1\n"); 29447c6ae99SBarry Smith } 2951321219cSEthan Coon if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)){ 29647c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\ 29747c6ae99SBarry Smith by 2*stencil_width + 1\n"); 29847c6ae99SBarry Smith } 29947c6ae99SBarry Smith 300aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 301aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 30247c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 30347c6ae99SBarry Smith 30447c6ae99SBarry Smith /* create the coloring */ 30547c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 30647c6ae99SBarry Smith if (!dd->localcoloring) { 30747c6ae99SBarry Smith ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 30847c6ae99SBarry Smith ii = 0; 30947c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 31047c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 31147c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 31247c6ae99SBarry Smith for (l=0; l<nc; l++) { 31347c6ae99SBarry Smith colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col)); 31447c6ae99SBarry Smith } 31547c6ae99SBarry Smith } 31647c6ae99SBarry Smith } 31747c6ae99SBarry Smith } 31847c6ae99SBarry Smith ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1)); 31947c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr); 32047c6ae99SBarry Smith } 32147c6ae99SBarry Smith *coloring = dd->localcoloring; 32247c6ae99SBarry Smith } else if (ctype == IS_COLORING_GHOSTED) { 32347c6ae99SBarry Smith if (!dd->ghostedcoloring) { 32447c6ae99SBarry Smith ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 32547c6ae99SBarry Smith ii = 0; 32647c6ae99SBarry Smith for (k=gzs; k<gzs+gnz; k++) { 32747c6ae99SBarry Smith for (j=gys; j<gys+gny; j++) { 32847c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 32947c6ae99SBarry Smith for (l=0; l<nc; l++) { 33047c6ae99SBarry Smith /* the complicated stuff is to handle periodic boundaries */ 33147c6ae99SBarry Smith colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col)); 33247c6ae99SBarry Smith } 33347c6ae99SBarry Smith } 33447c6ae99SBarry Smith } 33547c6ae99SBarry Smith } 33647c6ae99SBarry Smith ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1)); 33747c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr); 33847c6ae99SBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr); 33947c6ae99SBarry Smith } 34047c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 34147c6ae99SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 34247c6ae99SBarry Smith ierr = ISColoringReference(*coloring);CHKERRQ(ierr); 34347c6ae99SBarry Smith PetscFunctionReturn(0); 34447c6ae99SBarry Smith } 34547c6ae99SBarry Smith 34647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 34747c6ae99SBarry Smith 34847c6ae99SBarry Smith #undef __FUNCT__ 349e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_1d_MPIAIJ" 350e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 35147c6ae99SBarry Smith { 35247c6ae99SBarry Smith PetscErrorCode ierr; 35347c6ae99SBarry Smith PetscInt xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col; 35447c6ae99SBarry Smith PetscInt ncolors; 35547c6ae99SBarry Smith MPI_Comm comm; 3561321219cSEthan Coon DMDABoundaryType bx; 35747c6ae99SBarry Smith ISColoringValue *colors; 35847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 35947c6ae99SBarry Smith 36047c6ae99SBarry Smith PetscFunctionBegin; 36147c6ae99SBarry Smith /* 36247c6ae99SBarry Smith nc - number of components per grid point 36347c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 36447c6ae99SBarry Smith 36547c6ae99SBarry Smith */ 3661321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr); 36747c6ae99SBarry Smith col = 2*s + 1; 36847c6ae99SBarry Smith 3691321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) { 37031e6f798SBarry Smith SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points %d is divisible\n\ 37131e6f798SBarry Smith by 2*stencil_width + 1 %d\n",(int)m,(int)col); 37247c6ae99SBarry Smith } 37347c6ae99SBarry Smith 374aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr); 375aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr); 37647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 37747c6ae99SBarry Smith 37847c6ae99SBarry Smith /* create the coloring */ 37947c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 38047c6ae99SBarry Smith if (!dd->localcoloring) { 38147c6ae99SBarry Smith ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 38247c6ae99SBarry Smith i1 = 0; 38347c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 38447c6ae99SBarry Smith for (l=0; l<nc; l++) { 38547c6ae99SBarry Smith colors[i1++] = l + nc*(i % col); 38647c6ae99SBarry Smith } 38747c6ae99SBarry Smith } 38847c6ae99SBarry Smith ncolors = nc + nc*(col-1); 38947c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr); 39047c6ae99SBarry Smith } 39147c6ae99SBarry Smith *coloring = dd->localcoloring; 39247c6ae99SBarry Smith } else if (ctype == IS_COLORING_GHOSTED) { 39347c6ae99SBarry Smith if (!dd->ghostedcoloring) { 39447c6ae99SBarry Smith ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 39547c6ae99SBarry Smith i1 = 0; 39647c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 39747c6ae99SBarry Smith for (l=0; l<nc; l++) { 39847c6ae99SBarry Smith /* the complicated stuff is to handle periodic boundaries */ 39947c6ae99SBarry Smith colors[i1++] = l + nc*(SetInRange(i,m) % col); 40047c6ae99SBarry Smith } 40147c6ae99SBarry Smith } 40247c6ae99SBarry Smith ncolors = nc + nc*(col-1); 40347c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr); 40447c6ae99SBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr); 40547c6ae99SBarry Smith } 40647c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 40747c6ae99SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 40847c6ae99SBarry Smith ierr = ISColoringReference(*coloring);CHKERRQ(ierr); 40947c6ae99SBarry Smith PetscFunctionReturn(0); 41047c6ae99SBarry Smith } 41147c6ae99SBarry Smith 41247c6ae99SBarry Smith #undef __FUNCT__ 413e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_5pt_MPIAIJ" 414e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 41547c6ae99SBarry Smith { 41647c6ae99SBarry Smith PetscErrorCode ierr; 41747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc; 41847c6ae99SBarry Smith PetscInt ncolors; 41947c6ae99SBarry Smith MPI_Comm comm; 4201321219cSEthan Coon DMDABoundaryType bx,by; 42147c6ae99SBarry Smith ISColoringValue *colors; 42247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 42347c6ae99SBarry Smith 42447c6ae99SBarry Smith PetscFunctionBegin; 42547c6ae99SBarry Smith /* 42647c6ae99SBarry Smith nc - number of components per grid point 42747c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 42847c6ae99SBarry Smith 42947c6ae99SBarry Smith */ 4301321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr); 431aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 432aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 43347c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 43447c6ae99SBarry Smith 4351321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC && (m % 5)){ 43647c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\ 43747c6ae99SBarry Smith by 5\n"); 43847c6ae99SBarry Smith } 4391321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC && (n % 5)){ 44047c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\ 44147c6ae99SBarry Smith by 5\n"); 44247c6ae99SBarry Smith } 44347c6ae99SBarry Smith 44447c6ae99SBarry Smith /* create the coloring */ 44547c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 44647c6ae99SBarry Smith if (!dd->localcoloring) { 44747c6ae99SBarry Smith ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 44847c6ae99SBarry Smith ii = 0; 44947c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 45047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 45147c6ae99SBarry Smith for (k=0; k<nc; k++) { 45247c6ae99SBarry Smith colors[ii++] = k + nc*((3*j+i) % 5); 45347c6ae99SBarry Smith } 45447c6ae99SBarry Smith } 45547c6ae99SBarry Smith } 45647c6ae99SBarry Smith ncolors = 5*nc; 45747c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr); 45847c6ae99SBarry Smith } 45947c6ae99SBarry Smith *coloring = dd->localcoloring; 46047c6ae99SBarry Smith } else if (ctype == IS_COLORING_GHOSTED) { 46147c6ae99SBarry Smith if (!dd->ghostedcoloring) { 46247c6ae99SBarry Smith ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 46347c6ae99SBarry Smith ii = 0; 46447c6ae99SBarry Smith for (j=gys; j<gys+gny; j++) { 46547c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 46647c6ae99SBarry Smith for (k=0; k<nc; k++) { 46747c6ae99SBarry Smith colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5); 46847c6ae99SBarry Smith } 46947c6ae99SBarry Smith } 47047c6ae99SBarry Smith } 47147c6ae99SBarry Smith ncolors = 5*nc; 47247c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr); 47347c6ae99SBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr); 47447c6ae99SBarry Smith } 47547c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 47647c6ae99SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 47747c6ae99SBarry Smith PetscFunctionReturn(0); 47847c6ae99SBarry Smith } 47947c6ae99SBarry Smith 48047c6ae99SBarry Smith /* =========================================================================== */ 481950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM,Mat); 482950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM,Mat); 483950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM,Mat); 484950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM,Mat); 485950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM,Mat); 486950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM,Mat); 487950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM,Mat); 488950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM,Mat); 489950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM,Mat); 49047c6ae99SBarry Smith 49147c6ae99SBarry Smith #undef __FUNCT__ 49295ee5b0eSBarry Smith #define __FUNCT__ "MatSetDM" 49347c6ae99SBarry Smith /*@ 49495ee5b0eSBarry Smith MatSetDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix 49547c6ae99SBarry Smith 49647c6ae99SBarry Smith Logically Collective on Mat 49747c6ae99SBarry Smith 49847c6ae99SBarry Smith Input Parameters: 49947c6ae99SBarry Smith + mat - the matrix 50047c6ae99SBarry Smith - da - the da 50147c6ae99SBarry Smith 50247c6ae99SBarry Smith Level: intermediate 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith @*/ 50595ee5b0eSBarry Smith PetscErrorCode MatSetDM(Mat mat,DM da) 50647c6ae99SBarry Smith { 50747c6ae99SBarry Smith PetscErrorCode ierr; 50847c6ae99SBarry Smith 50947c6ae99SBarry Smith PetscFunctionBegin; 51047c6ae99SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 51147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 51295ee5b0eSBarry Smith ierr = PetscTryMethod(mat,"MatSetDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr); 51347c6ae99SBarry Smith PetscFunctionReturn(0); 51447c6ae99SBarry Smith } 51547c6ae99SBarry Smith 51647c6ae99SBarry Smith EXTERN_C_BEGIN 51747c6ae99SBarry Smith #undef __FUNCT__ 51847c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA" 5197087cfbeSBarry Smith PetscErrorCode MatView_MPI_DA(Mat A,PetscViewer viewer) 52047c6ae99SBarry Smith { 5219a42bb27SBarry Smith DM da; 52247c6ae99SBarry Smith PetscErrorCode ierr; 52347c6ae99SBarry Smith const char *prefix; 52447c6ae99SBarry Smith Mat Anatural; 52547c6ae99SBarry Smith AO ao; 52647c6ae99SBarry Smith PetscInt rstart,rend,*petsc,i; 52747c6ae99SBarry Smith IS is; 52847c6ae99SBarry Smith MPI_Comm comm; 52947c6ae99SBarry Smith 53047c6ae99SBarry Smith PetscFunctionBegin; 53147c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 5323c0c59f3SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"DM",(PetscObject*)&da);CHKERRQ(ierr); 533aa219208SBarry Smith if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA"); 53447c6ae99SBarry Smith 535aa219208SBarry Smith ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr); 53647c6ae99SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 53747c6ae99SBarry Smith ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr); 53847c6ae99SBarry Smith for (i=rstart; i<rend; i++) petsc[i-rstart] = i; 53947c6ae99SBarry Smith ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr); 54047c6ae99SBarry Smith ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 54147c6ae99SBarry Smith 54247c6ae99SBarry Smith /* call viewer on natural ordering */ 54347c6ae99SBarry Smith ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr); 544fcfd50ebSBarry Smith ierr = ISDestroy(&is);CHKERRQ(ierr); 54547c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr); 54647c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr); 54747c6ae99SBarry Smith ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr); 54847c6ae99SBarry Smith ierr = MatView(Anatural,viewer);CHKERRQ(ierr); 549fcfd50ebSBarry Smith ierr = MatDestroy(&Anatural);CHKERRQ(ierr); 55047c6ae99SBarry Smith PetscFunctionReturn(0); 55147c6ae99SBarry Smith } 55247c6ae99SBarry Smith EXTERN_C_END 55347c6ae99SBarry Smith 55447c6ae99SBarry Smith EXTERN_C_BEGIN 55547c6ae99SBarry Smith #undef __FUNCT__ 55647c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA" 5577087cfbeSBarry Smith PetscErrorCode MatLoad_MPI_DA(Mat A,PetscViewer viewer) 55847c6ae99SBarry Smith { 5599a42bb27SBarry Smith DM da; 56047c6ae99SBarry Smith PetscErrorCode ierr; 56147c6ae99SBarry Smith Mat Anatural,Aapp; 56247c6ae99SBarry Smith AO ao; 56347c6ae99SBarry Smith PetscInt rstart,rend,*app,i; 56447c6ae99SBarry Smith IS is; 56547c6ae99SBarry Smith MPI_Comm comm; 56647c6ae99SBarry Smith 56747c6ae99SBarry Smith PetscFunctionBegin; 56847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 5693c0c59f3SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"DM",(PetscObject*)&da);CHKERRQ(ierr); 570aa219208SBarry Smith if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA"); 57147c6ae99SBarry Smith 57247c6ae99SBarry Smith /* Load the matrix in natural ordering */ 57347c6ae99SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr); 57447c6ae99SBarry Smith ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr); 57547c6ae99SBarry Smith ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 57647c6ae99SBarry Smith ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr); 57747c6ae99SBarry Smith 57847c6ae99SBarry Smith /* Map natural ordering to application ordering and create IS */ 579aa219208SBarry Smith ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr); 58047c6ae99SBarry Smith ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr); 58147c6ae99SBarry Smith ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr); 58247c6ae99SBarry Smith for (i=rstart; i<rend; i++) app[i-rstart] = i; 58347c6ae99SBarry Smith ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr); 58447c6ae99SBarry Smith ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 58547c6ae99SBarry Smith 58647c6ae99SBarry Smith /* Do permutation and replace header */ 58747c6ae99SBarry Smith ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr); 58847c6ae99SBarry Smith ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr); 589fcfd50ebSBarry Smith ierr = ISDestroy(&is);CHKERRQ(ierr); 590fcfd50ebSBarry Smith ierr = MatDestroy(&Anatural);CHKERRQ(ierr); 59147c6ae99SBarry Smith PetscFunctionReturn(0); 59247c6ae99SBarry Smith } 59347c6ae99SBarry Smith EXTERN_C_END 59447c6ae99SBarry Smith 59547c6ae99SBarry Smith #undef __FUNCT__ 596950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA" 597950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA(DM da, const MatType mtype,Mat *J) 59847c6ae99SBarry Smith { 59947c6ae99SBarry Smith PetscErrorCode ierr; 60047c6ae99SBarry Smith PetscInt dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P; 60147c6ae99SBarry Smith Mat A; 60247c6ae99SBarry Smith MPI_Comm comm; 60347c6ae99SBarry Smith const MatType Atype; 60447c6ae99SBarry Smith void (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL; 60547c6ae99SBarry Smith MatType ttype[256]; 60647c6ae99SBarry Smith PetscBool flg; 60747c6ae99SBarry Smith PetscMPIInt size; 60847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 60947c6ae99SBarry Smith 61047c6ae99SBarry Smith PetscFunctionBegin; 61147c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES 61247c6ae99SBarry Smith ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 61347c6ae99SBarry Smith #endif 6145da5aae0SJed Brown if (!mtype) mtype = MATAIJ; 61547c6ae99SBarry Smith ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr); 616aa219208SBarry Smith ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA options","Mat");CHKERRQ(ierr); 617dd85299cSBarry Smith ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr); 61847c6ae99SBarry Smith ierr = PetscOptionsEnd(); 61947c6ae99SBarry Smith 62047c6ae99SBarry Smith /* 62147c6ae99SBarry Smith m 62247c6ae99SBarry Smith ------------------------------------------------------ 62347c6ae99SBarry Smith | | 62447c6ae99SBarry Smith | | 62547c6ae99SBarry Smith | ---------------------- | 62647c6ae99SBarry Smith | | | | 62747c6ae99SBarry Smith n | ny | | | 62847c6ae99SBarry Smith | | | | 62947c6ae99SBarry Smith | .--------------------- | 63047c6ae99SBarry Smith | (xs,ys) nx | 63147c6ae99SBarry Smith | . | 63247c6ae99SBarry Smith | (gxs,gys) | 63347c6ae99SBarry Smith | | 63447c6ae99SBarry Smith ----------------------------------------------------- 63547c6ae99SBarry Smith */ 63647c6ae99SBarry Smith 63747c6ae99SBarry Smith /* 63847c6ae99SBarry Smith nc - number of components per grid point 63947c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 64047c6ae99SBarry Smith 64147c6ae99SBarry Smith */ 6421321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 643aa219208SBarry Smith ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr); 64447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 64547c6ae99SBarry Smith ierr = MatCreate(comm,&A);CHKERRQ(ierr); 64647c6ae99SBarry Smith ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr); 64747c6ae99SBarry Smith ierr = MatSetType(A,(const MatType)ttype);CHKERRQ(ierr); 64895ee5b0eSBarry Smith ierr = MatSetDM(A,da);CHKERRQ(ierr); 64947c6ae99SBarry Smith ierr = MatSetFromOptions(A);CHKERRQ(ierr); 65047c6ae99SBarry Smith ierr = MatGetType(A,&Atype);CHKERRQ(ierr); 65147c6ae99SBarry Smith /* 652aa219208SBarry Smith We do not provide a getmatrix function in the DMDA operations because 653aa219208SBarry Smith the basic DMDA does not know about matrices. We think of DMDA as being more 65447c6ae99SBarry Smith more low-level than matrices. This is kind of cheating but, cause sometimes 655aa219208SBarry Smith we think of DMDA has higher level than matrices. 65647c6ae99SBarry Smith 65747c6ae99SBarry Smith We could switch based on Atype (or mtype), but we do not since the 65847c6ae99SBarry Smith specialized setting routines depend only the particular preallocation 65947c6ae99SBarry Smith details of the matrix, not the type itself. 66047c6ae99SBarry Smith */ 66147c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 66247c6ae99SBarry Smith if (!aij) { 66347c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 66447c6ae99SBarry Smith } 66547c6ae99SBarry Smith if (!aij) { 66647c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr); 66747c6ae99SBarry Smith if (!baij) { 66847c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr); 66947c6ae99SBarry Smith } 67047c6ae99SBarry Smith if (!baij){ 67147c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr); 67247c6ae99SBarry Smith if (!sbaij) { 67347c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr); 67447c6ae99SBarry Smith } 67547c6ae99SBarry Smith } 67647c6ae99SBarry Smith } 67747c6ae99SBarry Smith if (aij) { 67847c6ae99SBarry Smith if (dim == 1) { 679950540a4SJed Brown ierr = DMCreateMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr); 68047c6ae99SBarry Smith } else if (dim == 2) { 68147c6ae99SBarry Smith if (dd->ofill) { 682950540a4SJed Brown ierr = DMCreateMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr); 68347c6ae99SBarry Smith } else { 684950540a4SJed Brown ierr = DMCreateMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr); 68547c6ae99SBarry Smith } 68647c6ae99SBarry Smith } else if (dim == 3) { 68747c6ae99SBarry Smith if (dd->ofill) { 688950540a4SJed Brown ierr = DMCreateMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr); 68947c6ae99SBarry Smith } else { 690950540a4SJed Brown ierr = DMCreateMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr); 69147c6ae99SBarry Smith } 69247c6ae99SBarry Smith } 69347c6ae99SBarry Smith } else if (baij) { 69447c6ae99SBarry Smith if (dim == 2) { 695950540a4SJed Brown ierr = DMCreateMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr); 69647c6ae99SBarry Smith } else if (dim == 3) { 697950540a4SJed Brown ierr = DMCreateMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr); 69847c6ae99SBarry Smith } else { 699b17742caSSean Farley SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \ 700b17742caSSean Farley "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim); 70147c6ae99SBarry Smith } 70247c6ae99SBarry Smith } else if (sbaij) { 70347c6ae99SBarry Smith if (dim == 2) { 704950540a4SJed Brown ierr = DMCreateMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr); 70547c6ae99SBarry Smith } else if (dim == 3) { 706950540a4SJed Brown ierr = DMCreateMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr); 70747c6ae99SBarry Smith } else { 708b17742caSSean Farley SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \ 709b17742caSSean Farley "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim); 71047c6ae99SBarry Smith } 711869776cdSLisandro Dalcin } else { 712869776cdSLisandro Dalcin ISLocalToGlobalMapping ltog,ltogb; 713869776cdSLisandro Dalcin ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 714869776cdSLisandro Dalcin ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 715869776cdSLisandro Dalcin ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr); 716869776cdSLisandro Dalcin ierr = MatSetLocalToGlobalMappingBlock(A,ltogb,ltogb);CHKERRQ(ierr); 71747c6ae99SBarry Smith } 718aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr); 71947c6ae99SBarry Smith ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr); 7203c0c59f3SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"DM",(PetscObject)da);CHKERRQ(ierr); 72147c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 72247c6ae99SBarry Smith if (size > 1) { 72347c6ae99SBarry Smith /* change viewer to display matrix in natural ordering */ 72447c6ae99SBarry Smith ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr); 72547c6ae99SBarry Smith ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr); 72647c6ae99SBarry Smith } 72747c6ae99SBarry Smith *J = A; 72847c6ae99SBarry Smith PetscFunctionReturn(0); 72947c6ae99SBarry Smith } 73047c6ae99SBarry Smith 73147c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 73247c6ae99SBarry Smith #undef __FUNCT__ 733950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ" 734950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM da,Mat J) 73547c6ae99SBarry Smith { 73647c6ae99SBarry Smith PetscErrorCode ierr; 73747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny,m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p; 73847c6ae99SBarry Smith PetscInt lstart,lend,pstart,pend,*dnz,*onz; 73947c6ae99SBarry Smith MPI_Comm comm; 74047c6ae99SBarry Smith PetscScalar *values; 7411321219cSEthan Coon DMDABoundaryType bx,by; 74247c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 743aa219208SBarry Smith DMDAStencilType st; 74447c6ae99SBarry Smith 74547c6ae99SBarry Smith PetscFunctionBegin; 74647c6ae99SBarry Smith /* 74747c6ae99SBarry Smith nc - number of components per grid point 74847c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 74947c6ae99SBarry Smith 75047c6ae99SBarry Smith */ 7511321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 75247c6ae99SBarry Smith col = 2*s + 1; 753aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 754aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 75547c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 75647c6ae99SBarry Smith 75747c6ae99SBarry Smith ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr); 7581411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 7591411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 76047c6ae99SBarry Smith 76147c6ae99SBarry Smith /* determine the matrix preallocation information */ 76247c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr); 76347c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 76447c6ae99SBarry Smith 7651321219cSEthan Coon pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 7661321219cSEthan Coon pend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 76747c6ae99SBarry Smith 76847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 76947c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 77047c6ae99SBarry Smith 7711321219cSEthan Coon lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 7721321219cSEthan Coon lend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 77347c6ae99SBarry Smith 77447c6ae99SBarry Smith cnt = 0; 77547c6ae99SBarry Smith for (k=0; k<nc; k++) { 77647c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 77747c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 778aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */ 77947c6ae99SBarry Smith cols[cnt++] = k + nc*(slot + gnx*l + p); 78047c6ae99SBarry Smith } 78147c6ae99SBarry Smith } 78247c6ae99SBarry Smith } 78347c6ae99SBarry Smith rows[k] = k + nc*(slot); 78447c6ae99SBarry Smith } 785784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 78647c6ae99SBarry Smith } 78747c6ae99SBarry Smith } 78847c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 78947c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 79047c6ae99SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 79147c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 79247c6ae99SBarry Smith 793784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 794784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 79547c6ae99SBarry Smith 79647c6ae99SBarry Smith /* 79747c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 79847c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 79947c6ae99SBarry Smith PETSc ordering. 80047c6ae99SBarry Smith */ 801fcfd50ebSBarry Smith if (!da->prealloc_only) { 80247c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 80347c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 80447c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 80547c6ae99SBarry Smith 8061321219cSEthan Coon pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 8071321219cSEthan Coon pend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 80847c6ae99SBarry Smith 80947c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 81047c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 81147c6ae99SBarry Smith 8121321219cSEthan Coon lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 8131321219cSEthan Coon lend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 81447c6ae99SBarry Smith 81547c6ae99SBarry Smith cnt = 0; 81647c6ae99SBarry Smith for (k=0; k<nc; k++) { 81747c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 81847c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 819aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */ 82047c6ae99SBarry Smith cols[cnt++] = k + nc*(slot + gnx*l + p); 82147c6ae99SBarry Smith } 82247c6ae99SBarry Smith } 82347c6ae99SBarry Smith } 82447c6ae99SBarry Smith rows[k] = k + nc*(slot); 82547c6ae99SBarry Smith } 82647c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 82747c6ae99SBarry Smith } 82847c6ae99SBarry Smith } 82947c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 83047c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 83147c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 83247c6ae99SBarry Smith } 83347c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 83447c6ae99SBarry Smith PetscFunctionReturn(0); 83547c6ae99SBarry Smith } 83647c6ae99SBarry Smith 83747c6ae99SBarry Smith #undef __FUNCT__ 838950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ_Fill" 839950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J) 84047c6ae99SBarry Smith { 84147c6ae99SBarry Smith PetscErrorCode ierr; 84247c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 84347c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,row,col,cnt,l,p; 84447c6ae99SBarry Smith PetscInt lstart,lend,pstart,pend,*dnz,*onz; 84547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 84647c6ae99SBarry Smith PetscInt ifill_col,*ofill = dd->ofill, *dfill = dd->dfill; 84747c6ae99SBarry Smith MPI_Comm comm; 84847c6ae99SBarry Smith PetscScalar *values; 8491321219cSEthan Coon DMDABoundaryType bx,by; 85047c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 851aa219208SBarry Smith DMDAStencilType st; 85247c6ae99SBarry Smith 85347c6ae99SBarry Smith PetscFunctionBegin; 85447c6ae99SBarry Smith /* 85547c6ae99SBarry Smith nc - number of components per grid point 85647c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 85747c6ae99SBarry Smith 85847c6ae99SBarry Smith */ 8591321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 86047c6ae99SBarry Smith col = 2*s + 1; 861aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 862aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 86347c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 86447c6ae99SBarry Smith 86547c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 8661411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 8671411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 86847c6ae99SBarry Smith 86947c6ae99SBarry Smith /* determine the matrix preallocation information */ 87047c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr); 87147c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 87247c6ae99SBarry Smith 8731321219cSEthan Coon pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 8741321219cSEthan Coon pend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 87547c6ae99SBarry Smith 87647c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 87747c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 87847c6ae99SBarry Smith 8791321219cSEthan Coon lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 8801321219cSEthan Coon lend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 88147c6ae99SBarry Smith 88247c6ae99SBarry Smith for (k=0; k<nc; k++) { 88347c6ae99SBarry Smith cnt = 0; 88447c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 88547c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 88647c6ae99SBarry Smith if (l || p) { 887aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star */ 88847c6ae99SBarry Smith for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) 88947c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p); 89047c6ae99SBarry Smith } 89147c6ae99SBarry Smith } else { 89247c6ae99SBarry Smith if (dfill) { 89347c6ae99SBarry Smith for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) 89447c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p); 89547c6ae99SBarry Smith } else { 89647c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 89747c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + gnx*l + p); 89847c6ae99SBarry Smith } 89947c6ae99SBarry Smith } 90047c6ae99SBarry Smith } 90147c6ae99SBarry Smith } 90247c6ae99SBarry Smith row = k + nc*(slot); 903784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 90447c6ae99SBarry Smith } 90547c6ae99SBarry Smith } 90647c6ae99SBarry Smith } 90747c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 90847c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 90947c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 910784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 911784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 91247c6ae99SBarry Smith 91347c6ae99SBarry Smith /* 91447c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 91547c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 91647c6ae99SBarry Smith PETSc ordering. 91747c6ae99SBarry Smith */ 918fcfd50ebSBarry Smith if (!da->prealloc_only) { 91947c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 92047c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 92147c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 92247c6ae99SBarry Smith 9231321219cSEthan Coon pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 9241321219cSEthan Coon pend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 92547c6ae99SBarry Smith 92647c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 92747c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 92847c6ae99SBarry Smith 9291321219cSEthan Coon lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 9301321219cSEthan Coon lend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 93147c6ae99SBarry Smith 93247c6ae99SBarry Smith for (k=0; k<nc; k++) { 93347c6ae99SBarry Smith cnt = 0; 93447c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 93547c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 93647c6ae99SBarry Smith if (l || p) { 937aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star */ 93847c6ae99SBarry Smith for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) 93947c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p); 94047c6ae99SBarry Smith } 94147c6ae99SBarry Smith } else { 94247c6ae99SBarry Smith if (dfill) { 94347c6ae99SBarry Smith for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) 94447c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p); 94547c6ae99SBarry Smith } else { 94647c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 94747c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + gnx*l + p); 94847c6ae99SBarry Smith } 94947c6ae99SBarry Smith } 95047c6ae99SBarry Smith } 95147c6ae99SBarry Smith } 95247c6ae99SBarry Smith row = k + nc*(slot); 95347c6ae99SBarry Smith ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 95447c6ae99SBarry Smith } 95547c6ae99SBarry Smith } 95647c6ae99SBarry Smith } 95747c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 95847c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 95947c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96047c6ae99SBarry Smith } 96147c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 96247c6ae99SBarry Smith PetscFunctionReturn(0); 96347c6ae99SBarry Smith } 96447c6ae99SBarry Smith 96547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 96647c6ae99SBarry Smith 96747c6ae99SBarry Smith #undef __FUNCT__ 968950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ" 969950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM da,Mat J) 97047c6ae99SBarry Smith { 97147c6ae99SBarry Smith PetscErrorCode ierr; 97247c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 97347c6ae99SBarry Smith PetscInt m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL; 97447c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 97547c6ae99SBarry Smith MPI_Comm comm; 97647c6ae99SBarry Smith PetscScalar *values; 9771321219cSEthan Coon DMDABoundaryType bx,by,bz; 97847c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 979aa219208SBarry Smith DMDAStencilType st; 98047c6ae99SBarry Smith 98147c6ae99SBarry Smith PetscFunctionBegin; 98247c6ae99SBarry Smith /* 98347c6ae99SBarry Smith nc - number of components per grid point 98447c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 98547c6ae99SBarry Smith 98647c6ae99SBarry Smith */ 9871321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 98847c6ae99SBarry Smith col = 2*s + 1; 98947c6ae99SBarry Smith 990aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 991aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 99247c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 99347c6ae99SBarry Smith 99447c6ae99SBarry Smith ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr); 9951411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 9961411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 99747c6ae99SBarry Smith 99847c6ae99SBarry Smith /* determine the matrix preallocation information */ 99947c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr); 100047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 10011321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 10021321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 100347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 10041321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 10051321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 100647c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 10071321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 10081321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 100947c6ae99SBarry Smith 101047c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 101147c6ae99SBarry Smith 101247c6ae99SBarry Smith cnt = 0; 101347c6ae99SBarry Smith for (l=0; l<nc; l++) { 101447c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 101547c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 101647c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1017aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 101847c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk); 101947c6ae99SBarry Smith } 102047c6ae99SBarry Smith } 102147c6ae99SBarry Smith } 102247c6ae99SBarry Smith } 102347c6ae99SBarry Smith rows[l] = l + nc*(slot); 102447c6ae99SBarry Smith } 1025784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 102647c6ae99SBarry Smith } 102747c6ae99SBarry Smith } 102847c6ae99SBarry Smith } 102947c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 103047c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 103147c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 103247c6ae99SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 1033784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1034784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 103547c6ae99SBarry Smith 103647c6ae99SBarry Smith /* 103747c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 103847c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 103947c6ae99SBarry Smith PETSc ordering. 104047c6ae99SBarry Smith */ 1041fcfd50ebSBarry Smith if (!da->prealloc_only) { 104247c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 104347c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 104447c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 10451321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 10461321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 104747c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 10481321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 10491321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 105047c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 10511321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 10521321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 105347c6ae99SBarry Smith 105447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 105547c6ae99SBarry Smith 105647c6ae99SBarry Smith cnt = 0; 105747c6ae99SBarry Smith for (l=0; l<nc; l++) { 105847c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 105947c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 106047c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1061aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 106247c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk); 106347c6ae99SBarry Smith } 106447c6ae99SBarry Smith } 106547c6ae99SBarry Smith } 106647c6ae99SBarry Smith } 106747c6ae99SBarry Smith rows[l] = l + nc*(slot); 106847c6ae99SBarry Smith } 106947c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 107047c6ae99SBarry Smith } 107147c6ae99SBarry Smith } 107247c6ae99SBarry Smith } 107347c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 107447c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 107547c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 107647c6ae99SBarry Smith } 107747c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 107847c6ae99SBarry Smith PetscFunctionReturn(0); 107947c6ae99SBarry Smith } 108047c6ae99SBarry Smith 108147c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 108247c6ae99SBarry Smith 108347c6ae99SBarry Smith #undef __FUNCT__ 1084950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ" 1085950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM da,Mat J) 108647c6ae99SBarry Smith { 108747c6ae99SBarry Smith PetscErrorCode ierr; 108847c6ae99SBarry Smith PetscInt xs,nx,i,i1,slot,gxs,gnx; 108947c6ae99SBarry Smith PetscInt m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l; 109047c6ae99SBarry Smith PetscInt istart,iend; 109147c6ae99SBarry Smith PetscScalar *values; 10921321219cSEthan Coon DMDABoundaryType bx; 109347c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 109447c6ae99SBarry Smith 109547c6ae99SBarry Smith PetscFunctionBegin; 109647c6ae99SBarry Smith /* 109747c6ae99SBarry Smith nc - number of components per grid point 109847c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 109947c6ae99SBarry Smith 110047c6ae99SBarry Smith */ 11011321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr); 110247c6ae99SBarry Smith col = 2*s + 1; 110347c6ae99SBarry Smith 1104aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr); 1105aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr); 110647c6ae99SBarry Smith 110747c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr); 110847c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr); 110947c6ae99SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 111047c6ae99SBarry Smith ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr); 111147c6ae99SBarry Smith 11121411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 11131411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 1114784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1115784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 111647c6ae99SBarry Smith 111747c6ae99SBarry Smith /* 111847c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 111947c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 112047c6ae99SBarry Smith PETSc ordering. 112147c6ae99SBarry Smith */ 1122fcfd50ebSBarry Smith if (!da->prealloc_only) { 112347c6ae99SBarry Smith ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 112447c6ae99SBarry Smith ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 112547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 112647c6ae99SBarry Smith istart = PetscMax(-s,gxs - i); 112747c6ae99SBarry Smith iend = PetscMin(s,gxs + gnx - i - 1); 112847c6ae99SBarry Smith slot = i - gxs; 112947c6ae99SBarry Smith 113047c6ae99SBarry Smith cnt = 0; 113147c6ae99SBarry Smith for (l=0; l<nc; l++) { 113247c6ae99SBarry Smith for (i1=istart; i1<iend+1; i1++) { 113347c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + i1); 113447c6ae99SBarry Smith } 113547c6ae99SBarry Smith rows[l] = l + nc*(slot); 113647c6ae99SBarry Smith } 113747c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 113847c6ae99SBarry Smith } 113947c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 114047c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 114147c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 114247c6ae99SBarry Smith } 114347c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 114447c6ae99SBarry Smith PetscFunctionReturn(0); 114547c6ae99SBarry Smith } 114647c6ae99SBarry Smith 114747c6ae99SBarry Smith #undef __FUNCT__ 1148950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIBAIJ" 1149950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM da,Mat J) 115047c6ae99SBarry Smith { 115147c6ae99SBarry Smith PetscErrorCode ierr; 115247c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 115347c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz; 115447c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,ii,jj; 115547c6ae99SBarry Smith MPI_Comm comm; 115647c6ae99SBarry Smith PetscScalar *values; 11571321219cSEthan Coon DMDABoundaryType bx,by; 1158aa219208SBarry Smith DMDAStencilType st; 115947c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 116047c6ae99SBarry Smith 116147c6ae99SBarry Smith PetscFunctionBegin; 116247c6ae99SBarry Smith /* 116347c6ae99SBarry Smith nc - number of components per grid point 116447c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 116547c6ae99SBarry Smith */ 11661321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 116747c6ae99SBarry Smith col = 2*s + 1; 116847c6ae99SBarry Smith 1169aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 1170aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 117147c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 117247c6ae99SBarry Smith 117347c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 117447c6ae99SBarry Smith 11751411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 11761411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 117747c6ae99SBarry Smith 117847c6ae99SBarry Smith /* determine the matrix preallocation information */ 117947c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr); 118047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 11811321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 11821321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 118347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 11841321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 11851321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 118647c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 118747c6ae99SBarry Smith 118847c6ae99SBarry Smith /* Find block columns in block row */ 118947c6ae99SBarry Smith cnt = 0; 119047c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 119147c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1192aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */ 119347c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 119447c6ae99SBarry Smith } 119547c6ae99SBarry Smith } 119647c6ae99SBarry Smith } 1197784ac674SJed Brown ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr); 119847c6ae99SBarry Smith } 119947c6ae99SBarry Smith } 120047c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 120147c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 120247c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 120347c6ae99SBarry Smith 1204784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1205784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 120647c6ae99SBarry Smith 120747c6ae99SBarry Smith /* 120847c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 120947c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 121047c6ae99SBarry Smith PETSc ordering. 121147c6ae99SBarry Smith */ 1212fcfd50ebSBarry Smith if (!da->prealloc_only) { 121347c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 121447c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 121547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 12161321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 12171321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 121847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 12191321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 12201321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 122147c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 122247c6ae99SBarry Smith cnt = 0; 122347c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 122447c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1225aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */ 122647c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 122747c6ae99SBarry Smith } 122847c6ae99SBarry Smith } 122947c6ae99SBarry Smith } 123047c6ae99SBarry Smith ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 123147c6ae99SBarry Smith } 123247c6ae99SBarry Smith } 123347c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 123447c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 123547c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 123647c6ae99SBarry Smith } 123747c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 123847c6ae99SBarry Smith PetscFunctionReturn(0); 123947c6ae99SBarry Smith } 124047c6ae99SBarry Smith 124147c6ae99SBarry Smith #undef __FUNCT__ 1242950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIBAIJ" 1243950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM da,Mat J) 124447c6ae99SBarry Smith { 124547c6ae99SBarry Smith PetscErrorCode ierr; 124647c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 124747c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz; 124847c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 124947c6ae99SBarry Smith MPI_Comm comm; 125047c6ae99SBarry Smith PetscScalar *values; 12511321219cSEthan Coon DMDABoundaryType bx,by,bz; 1252aa219208SBarry Smith DMDAStencilType st; 125347c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 125447c6ae99SBarry Smith 125547c6ae99SBarry Smith PetscFunctionBegin; 125647c6ae99SBarry Smith /* 125747c6ae99SBarry Smith nc - number of components per grid point 125847c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 125947c6ae99SBarry Smith 126047c6ae99SBarry Smith */ 12611321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 126247c6ae99SBarry Smith col = 2*s + 1; 126347c6ae99SBarry Smith 1264aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1265aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 126647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 126747c6ae99SBarry Smith 126847c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr); 126947c6ae99SBarry Smith 12701411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 12711411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 127247c6ae99SBarry Smith 127347c6ae99SBarry Smith /* determine the matrix preallocation information */ 127447c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr); 127547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 12761321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 12771321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 127847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 12791321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 12801321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 128147c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 12821321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 12831321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 128447c6ae99SBarry Smith 128547c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 128647c6ae99SBarry Smith 128747c6ae99SBarry Smith /* Find block columns in block row */ 128847c6ae99SBarry Smith cnt = 0; 128947c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 129047c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 129147c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1292aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 129347c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 129447c6ae99SBarry Smith } 129547c6ae99SBarry Smith } 129647c6ae99SBarry Smith } 129747c6ae99SBarry Smith } 1298784ac674SJed Brown ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr); 129947c6ae99SBarry Smith } 130047c6ae99SBarry Smith } 130147c6ae99SBarry Smith } 130247c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 130347c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 130447c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 130547c6ae99SBarry Smith 1306784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1307784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 130847c6ae99SBarry Smith 130947c6ae99SBarry Smith /* 131047c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 131147c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 131247c6ae99SBarry Smith PETSc ordering. 131347c6ae99SBarry Smith */ 1314fcfd50ebSBarry Smith if (!da->prealloc_only) { 131547c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 131647c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 131747c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 13181321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 13191321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 132047c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 13211321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 13221321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 132347c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 13241321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 13251321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 132647c6ae99SBarry Smith 132747c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 132847c6ae99SBarry Smith 132947c6ae99SBarry Smith cnt = 0; 133047c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 133147c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 133247c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1333aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 133447c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 133547c6ae99SBarry Smith } 133647c6ae99SBarry Smith } 133747c6ae99SBarry Smith } 133847c6ae99SBarry Smith } 133947c6ae99SBarry Smith ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 134047c6ae99SBarry Smith } 134147c6ae99SBarry Smith } 134247c6ae99SBarry Smith } 134347c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 134447c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 134547c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 134647c6ae99SBarry Smith } 134747c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 134847c6ae99SBarry Smith PetscFunctionReturn(0); 134947c6ae99SBarry Smith } 135047c6ae99SBarry Smith 135147c6ae99SBarry Smith #undef __FUNCT__ 135247c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular" 135347c6ae99SBarry Smith /* 135447c6ae99SBarry Smith This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to 135547c6ae99SBarry Smith identify in the local ordering with periodic domain. 135647c6ae99SBarry Smith */ 135747c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[]) 135847c6ae99SBarry Smith { 135947c6ae99SBarry Smith PetscErrorCode ierr; 136047c6ae99SBarry Smith PetscInt i,n; 136147c6ae99SBarry Smith 136247c6ae99SBarry Smith PetscFunctionBegin; 136347c6ae99SBarry Smith ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr); 136447c6ae99SBarry Smith ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr); 136547c6ae99SBarry Smith for (i=0,n=0; i<*cnt; i++) { 136647c6ae99SBarry Smith if (col[i] >= *row) col[n++] = col[i]; 136747c6ae99SBarry Smith } 136847c6ae99SBarry Smith *cnt = n; 136947c6ae99SBarry Smith PetscFunctionReturn(0); 137047c6ae99SBarry Smith } 137147c6ae99SBarry Smith 137247c6ae99SBarry Smith #undef __FUNCT__ 1373950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPISBAIJ" 1374950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM da,Mat J) 137547c6ae99SBarry Smith { 137647c6ae99SBarry Smith PetscErrorCode ierr; 137747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 137847c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz; 137947c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,ii,jj; 138047c6ae99SBarry Smith MPI_Comm comm; 138147c6ae99SBarry Smith PetscScalar *values; 13821321219cSEthan Coon DMDABoundaryType bx,by; 1383aa219208SBarry Smith DMDAStencilType st; 138447c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 138547c6ae99SBarry Smith 138647c6ae99SBarry Smith PetscFunctionBegin; 138747c6ae99SBarry Smith /* 138847c6ae99SBarry Smith nc - number of components per grid point 138947c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 139047c6ae99SBarry Smith */ 13911321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 139247c6ae99SBarry Smith col = 2*s + 1; 139347c6ae99SBarry Smith 1394aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 1395aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 139647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 139747c6ae99SBarry Smith 139847c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 139947c6ae99SBarry Smith 14001411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 14011411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 140247c6ae99SBarry Smith 140347c6ae99SBarry Smith /* determine the matrix preallocation information */ 140447c6ae99SBarry Smith ierr = MatPreallocateSymmetricInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr); 140547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 14061321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 14071321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 140847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 14091321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 14101321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 141147c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 141247c6ae99SBarry Smith 141347c6ae99SBarry Smith /* Find block columns in block row */ 141447c6ae99SBarry Smith cnt = 0; 141547c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 141647c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1417aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { 141847c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 141947c6ae99SBarry Smith } 142047c6ae99SBarry Smith } 142147c6ae99SBarry Smith } 142247c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 142347c6ae99SBarry Smith ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr); 142447c6ae99SBarry Smith } 142547c6ae99SBarry Smith } 142647c6ae99SBarry Smith ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 142747c6ae99SBarry Smith ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 142847c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 142947c6ae99SBarry Smith 1430784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1431784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 143247c6ae99SBarry Smith 143347c6ae99SBarry Smith /* 143447c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 143547c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 143647c6ae99SBarry Smith PETSc ordering. 143747c6ae99SBarry Smith */ 1438fcfd50ebSBarry Smith if (!da->prealloc_only) { 143947c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 144047c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 144147c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 14421321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 14431321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 144447c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 14451321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 14461321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 144747c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 144847c6ae99SBarry Smith 144947c6ae99SBarry Smith /* Find block columns in block row */ 145047c6ae99SBarry Smith cnt = 0; 145147c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 145247c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1453aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { 145447c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 145547c6ae99SBarry Smith } 145647c6ae99SBarry Smith } 145747c6ae99SBarry Smith } 145847c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 145947c6ae99SBarry Smith ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 146047c6ae99SBarry Smith } 146147c6ae99SBarry Smith } 146247c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 146347c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 146447c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 146547c6ae99SBarry Smith } 146647c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 146747c6ae99SBarry Smith PetscFunctionReturn(0); 146847c6ae99SBarry Smith } 146947c6ae99SBarry Smith 147047c6ae99SBarry Smith #undef __FUNCT__ 1471950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPISBAIJ" 1472950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM da,Mat J) 147347c6ae99SBarry Smith { 147447c6ae99SBarry Smith PetscErrorCode ierr; 147547c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 147647c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz; 147747c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 147847c6ae99SBarry Smith MPI_Comm comm; 147947c6ae99SBarry Smith PetscScalar *values; 14801321219cSEthan Coon DMDABoundaryType bx,by,bz; 1481aa219208SBarry Smith DMDAStencilType st; 148247c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 148347c6ae99SBarry Smith 148447c6ae99SBarry Smith PetscFunctionBegin; 148547c6ae99SBarry Smith /* 148647c6ae99SBarry Smith nc - number of components per grid point 148747c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 148847c6ae99SBarry Smith */ 14891321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 149047c6ae99SBarry Smith col = 2*s + 1; 149147c6ae99SBarry Smith 1492aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1493aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 149447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 149547c6ae99SBarry Smith 149647c6ae99SBarry Smith /* create the matrix */ 149747c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr); 149847c6ae99SBarry Smith 14991411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 15001411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 150147c6ae99SBarry Smith 150247c6ae99SBarry Smith /* determine the matrix preallocation information */ 150347c6ae99SBarry Smith ierr = MatPreallocateSymmetricInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr); 150447c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 15051321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 15061321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 150747c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 15081321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 15091321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 151047c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 15111321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 15121321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 151347c6ae99SBarry Smith 151447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 151547c6ae99SBarry Smith 151647c6ae99SBarry Smith /* Find block columns in block row */ 151747c6ae99SBarry Smith cnt = 0; 151847c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 151947c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 152047c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1521aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) { 152247c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 152347c6ae99SBarry Smith } 152447c6ae99SBarry Smith } 152547c6ae99SBarry Smith } 152647c6ae99SBarry Smith } 152747c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 152847c6ae99SBarry Smith ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr); 152947c6ae99SBarry Smith } 153047c6ae99SBarry Smith } 153147c6ae99SBarry Smith } 153247c6ae99SBarry Smith ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 153347c6ae99SBarry Smith ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 153447c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 153547c6ae99SBarry Smith 1536784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1537784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 153847c6ae99SBarry Smith 153947c6ae99SBarry Smith /* 154047c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 154147c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 154247c6ae99SBarry Smith PETSc ordering. 154347c6ae99SBarry Smith */ 1544fcfd50ebSBarry Smith if (!da->prealloc_only) { 154547c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 154647c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 154747c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 15481321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 15491321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 155047c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 15511321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 15521321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 155347c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 15541321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 15551321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 155647c6ae99SBarry Smith 155747c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 155847c6ae99SBarry Smith 155947c6ae99SBarry Smith cnt = 0; 156047c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 156147c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 156247c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1563aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) { 156447c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 156547c6ae99SBarry Smith } 156647c6ae99SBarry Smith } 156747c6ae99SBarry Smith } 156847c6ae99SBarry Smith } 156947c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 157047c6ae99SBarry Smith ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 157147c6ae99SBarry Smith } 157247c6ae99SBarry Smith } 157347c6ae99SBarry Smith } 157447c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 157547c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157647c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157747c6ae99SBarry Smith } 157847c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 157947c6ae99SBarry Smith PetscFunctionReturn(0); 158047c6ae99SBarry Smith } 158147c6ae99SBarry Smith 158247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 158347c6ae99SBarry Smith 158447c6ae99SBarry Smith #undef __FUNCT__ 1585950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ_Fill" 1586950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J) 158747c6ae99SBarry Smith { 158847c6ae99SBarry Smith PetscErrorCode ierr; 158947c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 159047c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz; 159147c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 159247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 159347c6ae99SBarry Smith PetscInt ifill_col,*dfill = dd->dfill,*ofill = dd->ofill; 159447c6ae99SBarry Smith MPI_Comm comm; 159547c6ae99SBarry Smith PetscScalar *values; 15961321219cSEthan Coon DMDABoundaryType bx,by,bz; 159747c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 1598aa219208SBarry Smith DMDAStencilType st; 159947c6ae99SBarry Smith 160047c6ae99SBarry Smith PetscFunctionBegin; 160147c6ae99SBarry Smith /* 160247c6ae99SBarry Smith nc - number of components per grid point 160347c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 160447c6ae99SBarry Smith 160547c6ae99SBarry Smith */ 16061321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 160747c6ae99SBarry Smith col = 2*s + 1; 16081321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){ 160947c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\ 161047c6ae99SBarry Smith by 2*stencil_width + 1\n"); 161147c6ae99SBarry Smith } 16121321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){ 161347c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\ 161447c6ae99SBarry Smith by 2*stencil_width + 1\n"); 161547c6ae99SBarry Smith } 16161321219cSEthan Coon if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)){ 161747c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\ 161847c6ae99SBarry Smith by 2*stencil_width + 1\n"); 161947c6ae99SBarry Smith } 162047c6ae99SBarry Smith 1621aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1622aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 162347c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 162447c6ae99SBarry Smith 162547c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 16261411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 16271411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 162847c6ae99SBarry Smith 162947c6ae99SBarry Smith /* determine the matrix preallocation information */ 163047c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr); 163147c6ae99SBarry Smith 163247c6ae99SBarry Smith 163347c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 16341321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 16351321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 163647c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 16371321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 16381321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 163947c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 16401321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 16411321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 164247c6ae99SBarry Smith 164347c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 164447c6ae99SBarry Smith 164547c6ae99SBarry Smith for (l=0; l<nc; l++) { 164647c6ae99SBarry Smith cnt = 0; 164747c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 164847c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 164947c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 165047c6ae99SBarry Smith if (ii || jj || kk) { 1651aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 165247c6ae99SBarry Smith for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) 165347c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 165447c6ae99SBarry Smith } 165547c6ae99SBarry Smith } else { 165647c6ae99SBarry Smith if (dfill) { 165747c6ae99SBarry Smith for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) 165847c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 165947c6ae99SBarry Smith } else { 166047c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 166147c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk); 166247c6ae99SBarry Smith } 166347c6ae99SBarry Smith } 166447c6ae99SBarry Smith } 166547c6ae99SBarry Smith } 166647c6ae99SBarry Smith } 166747c6ae99SBarry Smith row = l + nc*(slot); 1668784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 166947c6ae99SBarry Smith } 167047c6ae99SBarry Smith } 167147c6ae99SBarry Smith } 167247c6ae99SBarry Smith } 167347c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 167447c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 167547c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1676784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1677784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 167847c6ae99SBarry Smith 167947c6ae99SBarry Smith /* 168047c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 168147c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 168247c6ae99SBarry Smith PETSc ordering. 168347c6ae99SBarry Smith */ 1684fcfd50ebSBarry Smith if (!da->prealloc_only) { 168547c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 168647c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 168747c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 16881321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 16891321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 169047c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 16911321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 16921321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 169347c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 16941321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 16951321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 169647c6ae99SBarry Smith 169747c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 169847c6ae99SBarry Smith 169947c6ae99SBarry Smith for (l=0; l<nc; l++) { 170047c6ae99SBarry Smith cnt = 0; 170147c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 170247c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 170347c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 170447c6ae99SBarry Smith if (ii || jj || kk) { 1705aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 170647c6ae99SBarry Smith for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) 170747c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 170847c6ae99SBarry Smith } 170947c6ae99SBarry Smith } else { 171047c6ae99SBarry Smith if (dfill) { 171147c6ae99SBarry Smith for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) 171247c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 171347c6ae99SBarry Smith } else { 171447c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 171547c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk); 171647c6ae99SBarry Smith } 171747c6ae99SBarry Smith } 171847c6ae99SBarry Smith } 171947c6ae99SBarry Smith } 172047c6ae99SBarry Smith } 172147c6ae99SBarry Smith row = l + nc*(slot); 172247c6ae99SBarry Smith ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 172347c6ae99SBarry Smith } 172447c6ae99SBarry Smith } 172547c6ae99SBarry Smith } 172647c6ae99SBarry Smith } 172747c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 172847c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 172947c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 173047c6ae99SBarry Smith } 173147c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 173247c6ae99SBarry Smith PetscFunctionReturn(0); 173347c6ae99SBarry Smith } 1734