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 609573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring *); 709573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring *); 809573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring *); 909573ac7SBarry Smith extern PetscErrorCode DMGetColoring_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 5694013140SBarry Smith of the matrix returned by DMGetMatrix(). 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 85aa219208SBarry Smith .seealso DMGetMatrix(), DMDASetGetMatrix(), DMDASetMatPreallocateOnly() 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__ 10194013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA" 1027087cfbeSBarry Smith PetscErrorCode DMGetColoring_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) { 16894013140SBarry Smith ierr = DMGetColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 16947c6ae99SBarry Smith } else if (dim == 2) { 17094013140SBarry Smith ierr = DMGetColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 17147c6ae99SBarry Smith } else if (dim == 3) { 17294013140SBarry Smith ierr = DMGetColoring_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__ 18794013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_2d_MPIAIJ" 18894013140SBarry Smith PetscErrorCode DMGetColoring_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) { 21394013140SBarry Smith ierr = DMGetColoring_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__ 26794013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_3d_MPIAIJ" 26894013140SBarry Smith PetscErrorCode DMGetColoring_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__ 34994013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_1d_MPIAIJ" 35094013140SBarry Smith PetscErrorCode DMGetColoring_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__ 41394013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_2d_5pt_MPIAIJ" 41494013140SBarry Smith PetscErrorCode DMGetColoring_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 /* =========================================================================== */ 48109573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM,Mat); 48209573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM,Mat); 48309573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM,Mat); 48409573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM,Mat); 48509573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM,Mat); 48609573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM,Mat); 48709573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM,Mat); 48809573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM,Mat); 48909573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_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; 529*74388724SJed Brown PetscViewerFormat format; 53047c6ae99SBarry Smith 53147c6ae99SBarry Smith PetscFunctionBegin; 532*74388724SJed Brown /* Check whether we are just printing info, in which case MatView() already viewed everything we wanted to view */ 533*74388724SJed Brown ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 534*74388724SJed Brown if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0); 535*74388724SJed Brown 53647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 5373c0c59f3SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"DM",(PetscObject*)&da);CHKERRQ(ierr); 538aa219208SBarry Smith if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA"); 53947c6ae99SBarry Smith 540aa219208SBarry Smith ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr); 54147c6ae99SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 54247c6ae99SBarry Smith ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr); 54347c6ae99SBarry Smith for (i=rstart; i<rend; i++) petsc[i-rstart] = i; 54447c6ae99SBarry Smith ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr); 54547c6ae99SBarry Smith ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 54647c6ae99SBarry Smith 54747c6ae99SBarry Smith /* call viewer on natural ordering */ 54847c6ae99SBarry Smith ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr); 549fcfd50ebSBarry Smith ierr = ISDestroy(&is);CHKERRQ(ierr); 55047c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr); 55147c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr); 55247c6ae99SBarry Smith ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr); 55347c6ae99SBarry Smith ierr = MatView(Anatural,viewer);CHKERRQ(ierr); 554fcfd50ebSBarry Smith ierr = MatDestroy(&Anatural);CHKERRQ(ierr); 55547c6ae99SBarry Smith PetscFunctionReturn(0); 55647c6ae99SBarry Smith } 55747c6ae99SBarry Smith EXTERN_C_END 55847c6ae99SBarry Smith 55947c6ae99SBarry Smith EXTERN_C_BEGIN 56047c6ae99SBarry Smith #undef __FUNCT__ 56147c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA" 5627087cfbeSBarry Smith PetscErrorCode MatLoad_MPI_DA(Mat A,PetscViewer viewer) 56347c6ae99SBarry Smith { 5649a42bb27SBarry Smith DM da; 56547c6ae99SBarry Smith PetscErrorCode ierr; 56647c6ae99SBarry Smith Mat Anatural,Aapp; 56747c6ae99SBarry Smith AO ao; 56847c6ae99SBarry Smith PetscInt rstart,rend,*app,i; 56947c6ae99SBarry Smith IS is; 57047c6ae99SBarry Smith MPI_Comm comm; 57147c6ae99SBarry Smith 57247c6ae99SBarry Smith PetscFunctionBegin; 57347c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 5743c0c59f3SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"DM",(PetscObject*)&da);CHKERRQ(ierr); 575aa219208SBarry Smith if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA"); 57647c6ae99SBarry Smith 57747c6ae99SBarry Smith /* Load the matrix in natural ordering */ 57847c6ae99SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr); 57947c6ae99SBarry Smith ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr); 58047c6ae99SBarry Smith ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 58147c6ae99SBarry Smith ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr); 58247c6ae99SBarry Smith 58347c6ae99SBarry Smith /* Map natural ordering to application ordering and create IS */ 584aa219208SBarry Smith ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr); 58547c6ae99SBarry Smith ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr); 58647c6ae99SBarry Smith ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr); 58747c6ae99SBarry Smith for (i=rstart; i<rend; i++) app[i-rstart] = i; 58847c6ae99SBarry Smith ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr); 58947c6ae99SBarry Smith ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 59047c6ae99SBarry Smith 59147c6ae99SBarry Smith /* Do permutation and replace header */ 59247c6ae99SBarry Smith ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr); 59347c6ae99SBarry Smith ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr); 594fcfd50ebSBarry Smith ierr = ISDestroy(&is);CHKERRQ(ierr); 595fcfd50ebSBarry Smith ierr = MatDestroy(&Anatural);CHKERRQ(ierr); 59647c6ae99SBarry Smith PetscFunctionReturn(0); 59747c6ae99SBarry Smith } 59847c6ae99SBarry Smith EXTERN_C_END 59947c6ae99SBarry Smith 60047c6ae99SBarry Smith #undef __FUNCT__ 60194013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA" 6027087cfbeSBarry Smith PetscErrorCode DMGetMatrix_DA(DM da, const MatType mtype,Mat *J) 60347c6ae99SBarry Smith { 60447c6ae99SBarry Smith PetscErrorCode ierr; 60547c6ae99SBarry Smith PetscInt dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P; 60647c6ae99SBarry Smith Mat A; 60747c6ae99SBarry Smith MPI_Comm comm; 60847c6ae99SBarry Smith const MatType Atype; 60947c6ae99SBarry Smith void (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL; 61047c6ae99SBarry Smith MatType ttype[256]; 61147c6ae99SBarry Smith PetscBool flg; 61247c6ae99SBarry Smith PetscMPIInt size; 61347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 61447c6ae99SBarry Smith 61547c6ae99SBarry Smith PetscFunctionBegin; 61647c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES 61747c6ae99SBarry Smith ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 61847c6ae99SBarry Smith #endif 6195da5aae0SJed Brown if (!mtype) mtype = MATAIJ; 62047c6ae99SBarry Smith ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr); 621aa219208SBarry Smith ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA options","Mat");CHKERRQ(ierr); 62247c6ae99SBarry Smith ierr = PetscOptionsList("-da_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr); 62347c6ae99SBarry Smith ierr = PetscOptionsEnd(); 62447c6ae99SBarry Smith 62547c6ae99SBarry Smith /* 62647c6ae99SBarry Smith m 62747c6ae99SBarry Smith ------------------------------------------------------ 62847c6ae99SBarry Smith | | 62947c6ae99SBarry Smith | | 63047c6ae99SBarry Smith | ---------------------- | 63147c6ae99SBarry Smith | | | | 63247c6ae99SBarry Smith n | ny | | | 63347c6ae99SBarry Smith | | | | 63447c6ae99SBarry Smith | .--------------------- | 63547c6ae99SBarry Smith | (xs,ys) nx | 63647c6ae99SBarry Smith | . | 63747c6ae99SBarry Smith | (gxs,gys) | 63847c6ae99SBarry Smith | | 63947c6ae99SBarry Smith ----------------------------------------------------- 64047c6ae99SBarry Smith */ 64147c6ae99SBarry Smith 64247c6ae99SBarry Smith /* 64347c6ae99SBarry Smith nc - number of components per grid point 64447c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 64547c6ae99SBarry Smith 64647c6ae99SBarry Smith */ 6471321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 648aa219208SBarry Smith ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr); 64947c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 65047c6ae99SBarry Smith ierr = MatCreate(comm,&A);CHKERRQ(ierr); 65147c6ae99SBarry Smith ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr); 65247c6ae99SBarry Smith ierr = MatSetType(A,(const MatType)ttype);CHKERRQ(ierr); 65395ee5b0eSBarry Smith ierr = MatSetDM(A,da);CHKERRQ(ierr); 65447c6ae99SBarry Smith ierr = MatSetFromOptions(A);CHKERRQ(ierr); 65547c6ae99SBarry Smith ierr = MatGetType(A,&Atype);CHKERRQ(ierr); 65647c6ae99SBarry Smith /* 657aa219208SBarry Smith We do not provide a getmatrix function in the DMDA operations because 658aa219208SBarry Smith the basic DMDA does not know about matrices. We think of DMDA as being more 65947c6ae99SBarry Smith more low-level than matrices. This is kind of cheating but, cause sometimes 660aa219208SBarry Smith we think of DMDA has higher level than matrices. 66147c6ae99SBarry Smith 66247c6ae99SBarry Smith We could switch based on Atype (or mtype), but we do not since the 66347c6ae99SBarry Smith specialized setting routines depend only the particular preallocation 66447c6ae99SBarry Smith details of the matrix, not the type itself. 66547c6ae99SBarry Smith */ 66647c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 66747c6ae99SBarry Smith if (!aij) { 66847c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 66947c6ae99SBarry Smith } 67047c6ae99SBarry Smith if (!aij) { 67147c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr); 67247c6ae99SBarry Smith if (!baij) { 67347c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr); 67447c6ae99SBarry Smith } 67547c6ae99SBarry Smith if (!baij){ 67647c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr); 67747c6ae99SBarry Smith if (!sbaij) { 67847c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr); 67947c6ae99SBarry Smith } 68047c6ae99SBarry Smith } 68147c6ae99SBarry Smith } 68247c6ae99SBarry Smith if (aij) { 68347c6ae99SBarry Smith if (dim == 1) { 68494013140SBarry Smith ierr = DMGetMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr); 68547c6ae99SBarry Smith } else if (dim == 2) { 68647c6ae99SBarry Smith if (dd->ofill) { 68794013140SBarry Smith ierr = DMGetMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr); 68847c6ae99SBarry Smith } else { 68994013140SBarry Smith ierr = DMGetMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr); 69047c6ae99SBarry Smith } 69147c6ae99SBarry Smith } else if (dim == 3) { 69247c6ae99SBarry Smith if (dd->ofill) { 69394013140SBarry Smith ierr = DMGetMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr); 69447c6ae99SBarry Smith } else { 69594013140SBarry Smith ierr = DMGetMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr); 69647c6ae99SBarry Smith } 69747c6ae99SBarry Smith } 69847c6ae99SBarry Smith } else if (baij) { 69947c6ae99SBarry Smith if (dim == 2) { 70094013140SBarry Smith ierr = DMGetMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr); 70147c6ae99SBarry Smith } else if (dim == 3) { 70294013140SBarry Smith ierr = DMGetMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr); 70347c6ae99SBarry Smith } else { 704b17742caSSean Farley SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \ 705b17742caSSean Farley "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim); 70647c6ae99SBarry Smith } 70747c6ae99SBarry Smith } else if (sbaij) { 70847c6ae99SBarry Smith if (dim == 2) { 70994013140SBarry Smith ierr = DMGetMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr); 71047c6ae99SBarry Smith } else if (dim == 3) { 71194013140SBarry Smith ierr = DMGetMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr); 71247c6ae99SBarry Smith } else { 713b17742caSSean Farley SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \ 714b17742caSSean Farley "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim); 71547c6ae99SBarry Smith } 716869776cdSLisandro Dalcin } else { 717869776cdSLisandro Dalcin ISLocalToGlobalMapping ltog,ltogb; 718869776cdSLisandro Dalcin ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 719869776cdSLisandro Dalcin ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 720869776cdSLisandro Dalcin ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr); 721869776cdSLisandro Dalcin ierr = MatSetLocalToGlobalMappingBlock(A,ltogb,ltogb);CHKERRQ(ierr); 72247c6ae99SBarry Smith } 723aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr); 72447c6ae99SBarry Smith ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr); 7253c0c59f3SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"DM",(PetscObject)da);CHKERRQ(ierr); 72647c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 72747c6ae99SBarry Smith if (size > 1) { 72847c6ae99SBarry Smith /* change viewer to display matrix in natural ordering */ 72947c6ae99SBarry Smith ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr); 73047c6ae99SBarry Smith ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr); 73147c6ae99SBarry Smith } 73247c6ae99SBarry Smith *J = A; 73347c6ae99SBarry Smith PetscFunctionReturn(0); 73447c6ae99SBarry Smith } 73547c6ae99SBarry Smith 73647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 73747c6ae99SBarry Smith #undef __FUNCT__ 73894013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIAIJ" 73994013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM da,Mat J) 74047c6ae99SBarry Smith { 74147c6ae99SBarry Smith PetscErrorCode ierr; 74247c6ae99SBarry 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; 74347c6ae99SBarry Smith PetscInt lstart,lend,pstart,pend,*dnz,*onz; 74447c6ae99SBarry Smith MPI_Comm comm; 74547c6ae99SBarry Smith PetscScalar *values; 7461321219cSEthan Coon DMDABoundaryType bx,by; 74747c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 748aa219208SBarry Smith DMDAStencilType st; 74947c6ae99SBarry Smith 75047c6ae99SBarry Smith PetscFunctionBegin; 75147c6ae99SBarry Smith /* 75247c6ae99SBarry Smith nc - number of components per grid point 75347c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 75447c6ae99SBarry Smith 75547c6ae99SBarry Smith */ 7561321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 75747c6ae99SBarry Smith col = 2*s + 1; 758aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 759aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 76047c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 76147c6ae99SBarry Smith 76247c6ae99SBarry Smith ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr); 7631411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 7641411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 76547c6ae99SBarry Smith 76647c6ae99SBarry Smith /* determine the matrix preallocation information */ 76747c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr); 76847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 76947c6ae99SBarry Smith 7701321219cSEthan Coon pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 7711321219cSEthan Coon pend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 77247c6ae99SBarry Smith 77347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 77447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 77547c6ae99SBarry Smith 7761321219cSEthan Coon lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 7771321219cSEthan Coon lend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 77847c6ae99SBarry Smith 77947c6ae99SBarry Smith cnt = 0; 78047c6ae99SBarry Smith for (k=0; k<nc; k++) { 78147c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 78247c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 783aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */ 78447c6ae99SBarry Smith cols[cnt++] = k + nc*(slot + gnx*l + p); 78547c6ae99SBarry Smith } 78647c6ae99SBarry Smith } 78747c6ae99SBarry Smith } 78847c6ae99SBarry Smith rows[k] = k + nc*(slot); 78947c6ae99SBarry Smith } 790784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 79147c6ae99SBarry Smith } 79247c6ae99SBarry Smith } 79347c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 79447c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 79547c6ae99SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 79647c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 79747c6ae99SBarry Smith 798784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 799784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 80047c6ae99SBarry Smith 80147c6ae99SBarry Smith /* 80247c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 80347c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 80447c6ae99SBarry Smith PETSc ordering. 80547c6ae99SBarry Smith */ 806fcfd50ebSBarry Smith if (!da->prealloc_only) { 80747c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 80847c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 80947c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 81047c6ae99SBarry Smith 8111321219cSEthan Coon pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 8121321219cSEthan Coon pend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 81347c6ae99SBarry Smith 81447c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 81547c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 81647c6ae99SBarry Smith 8171321219cSEthan Coon lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 8181321219cSEthan Coon lend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 81947c6ae99SBarry Smith 82047c6ae99SBarry Smith cnt = 0; 82147c6ae99SBarry Smith for (k=0; k<nc; k++) { 82247c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 82347c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 824aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */ 82547c6ae99SBarry Smith cols[cnt++] = k + nc*(slot + gnx*l + p); 82647c6ae99SBarry Smith } 82747c6ae99SBarry Smith } 82847c6ae99SBarry Smith } 82947c6ae99SBarry Smith rows[k] = k + nc*(slot); 83047c6ae99SBarry Smith } 83147c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 83247c6ae99SBarry Smith } 83347c6ae99SBarry Smith } 83447c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 83547c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 83647c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 83747c6ae99SBarry Smith } 83847c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 83947c6ae99SBarry Smith PetscFunctionReturn(0); 84047c6ae99SBarry Smith } 84147c6ae99SBarry Smith 84247c6ae99SBarry Smith #undef __FUNCT__ 84394013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIAIJ_Fill" 84494013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J) 84547c6ae99SBarry Smith { 84647c6ae99SBarry Smith PetscErrorCode ierr; 84747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 84847c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,row,col,cnt,l,p; 84947c6ae99SBarry Smith PetscInt lstart,lend,pstart,pend,*dnz,*onz; 85047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 85147c6ae99SBarry Smith PetscInt ifill_col,*ofill = dd->ofill, *dfill = dd->dfill; 85247c6ae99SBarry Smith MPI_Comm comm; 85347c6ae99SBarry Smith PetscScalar *values; 8541321219cSEthan Coon DMDABoundaryType bx,by; 85547c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 856aa219208SBarry Smith DMDAStencilType st; 85747c6ae99SBarry Smith 85847c6ae99SBarry Smith PetscFunctionBegin; 85947c6ae99SBarry Smith /* 86047c6ae99SBarry Smith nc - number of components per grid point 86147c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 86247c6ae99SBarry Smith 86347c6ae99SBarry Smith */ 8641321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 86547c6ae99SBarry Smith col = 2*s + 1; 866aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 867aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 86847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 86947c6ae99SBarry Smith 87047c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 8711411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 8721411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 87347c6ae99SBarry Smith 87447c6ae99SBarry Smith /* determine the matrix preallocation information */ 87547c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr); 87647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 87747c6ae99SBarry Smith 8781321219cSEthan Coon pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 8791321219cSEthan Coon pend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 88047c6ae99SBarry Smith 88147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 88247c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 88347c6ae99SBarry Smith 8841321219cSEthan Coon lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 8851321219cSEthan Coon lend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 88647c6ae99SBarry Smith 88747c6ae99SBarry Smith for (k=0; k<nc; k++) { 88847c6ae99SBarry Smith cnt = 0; 88947c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 89047c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 89147c6ae99SBarry Smith if (l || p) { 892aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star */ 89347c6ae99SBarry Smith for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) 89447c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p); 89547c6ae99SBarry Smith } 89647c6ae99SBarry Smith } else { 89747c6ae99SBarry Smith if (dfill) { 89847c6ae99SBarry Smith for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) 89947c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p); 90047c6ae99SBarry Smith } else { 90147c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 90247c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + gnx*l + p); 90347c6ae99SBarry Smith } 90447c6ae99SBarry Smith } 90547c6ae99SBarry Smith } 90647c6ae99SBarry Smith } 90747c6ae99SBarry Smith row = k + nc*(slot); 908784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 90947c6ae99SBarry Smith } 91047c6ae99SBarry Smith } 91147c6ae99SBarry Smith } 91247c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 91347c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 91447c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 915784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 916784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 91747c6ae99SBarry Smith 91847c6ae99SBarry Smith /* 91947c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 92047c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 92147c6ae99SBarry Smith PETSc ordering. 92247c6ae99SBarry Smith */ 923fcfd50ebSBarry Smith if (!da->prealloc_only) { 92447c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 92547c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 92647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 92747c6ae99SBarry Smith 9281321219cSEthan Coon pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 9291321219cSEthan Coon pend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 93047c6ae99SBarry Smith 93147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 93247c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 93347c6ae99SBarry Smith 9341321219cSEthan Coon lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 9351321219cSEthan Coon lend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 93647c6ae99SBarry Smith 93747c6ae99SBarry Smith for (k=0; k<nc; k++) { 93847c6ae99SBarry Smith cnt = 0; 93947c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 94047c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 94147c6ae99SBarry Smith if (l || p) { 942aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star */ 94347c6ae99SBarry Smith for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) 94447c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p); 94547c6ae99SBarry Smith } 94647c6ae99SBarry Smith } else { 94747c6ae99SBarry Smith if (dfill) { 94847c6ae99SBarry Smith for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) 94947c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p); 95047c6ae99SBarry Smith } else { 95147c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 95247c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + gnx*l + p); 95347c6ae99SBarry Smith } 95447c6ae99SBarry Smith } 95547c6ae99SBarry Smith } 95647c6ae99SBarry Smith } 95747c6ae99SBarry Smith row = k + nc*(slot); 95847c6ae99SBarry Smith ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 95947c6ae99SBarry Smith } 96047c6ae99SBarry Smith } 96147c6ae99SBarry Smith } 96247c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 96347c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96447c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96547c6ae99SBarry Smith } 96647c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 96747c6ae99SBarry Smith PetscFunctionReturn(0); 96847c6ae99SBarry Smith } 96947c6ae99SBarry Smith 97047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 97147c6ae99SBarry Smith 97247c6ae99SBarry Smith #undef __FUNCT__ 97394013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIAIJ" 97494013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM da,Mat J) 97547c6ae99SBarry Smith { 97647c6ae99SBarry Smith PetscErrorCode ierr; 97747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 97847c6ae99SBarry Smith PetscInt m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL; 97947c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 98047c6ae99SBarry Smith MPI_Comm comm; 98147c6ae99SBarry Smith PetscScalar *values; 9821321219cSEthan Coon DMDABoundaryType bx,by,bz; 98347c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 984aa219208SBarry Smith DMDAStencilType st; 98547c6ae99SBarry Smith 98647c6ae99SBarry Smith PetscFunctionBegin; 98747c6ae99SBarry Smith /* 98847c6ae99SBarry Smith nc - number of components per grid point 98947c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 99047c6ae99SBarry Smith 99147c6ae99SBarry Smith */ 9921321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 99347c6ae99SBarry Smith col = 2*s + 1; 99447c6ae99SBarry Smith 995aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 996aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 99747c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 99847c6ae99SBarry Smith 99947c6ae99SBarry Smith ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr); 10001411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 10011411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 100247c6ae99SBarry Smith 100347c6ae99SBarry Smith /* determine the matrix preallocation information */ 100447c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr); 100547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 10061321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 10071321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 100847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 10091321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 10101321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 101147c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 10121321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 10131321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 101447c6ae99SBarry Smith 101547c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 101647c6ae99SBarry Smith 101747c6ae99SBarry Smith cnt = 0; 101847c6ae99SBarry Smith for (l=0; l<nc; l++) { 101947c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 102047c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 102147c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1022aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 102347c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk); 102447c6ae99SBarry Smith } 102547c6ae99SBarry Smith } 102647c6ae99SBarry Smith } 102747c6ae99SBarry Smith } 102847c6ae99SBarry Smith rows[l] = l + nc*(slot); 102947c6ae99SBarry Smith } 1030784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 103147c6ae99SBarry Smith } 103247c6ae99SBarry Smith } 103347c6ae99SBarry Smith } 103447c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 103547c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 103647c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 103747c6ae99SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 1038784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1039784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 104047c6ae99SBarry Smith 104147c6ae99SBarry Smith /* 104247c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 104347c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 104447c6ae99SBarry Smith PETSc ordering. 104547c6ae99SBarry Smith */ 1046fcfd50ebSBarry Smith if (!da->prealloc_only) { 104747c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 104847c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 104947c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 10501321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 10511321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 105247c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 10531321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 10541321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 105547c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 10561321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 10571321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 105847c6ae99SBarry Smith 105947c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 106047c6ae99SBarry Smith 106147c6ae99SBarry Smith cnt = 0; 106247c6ae99SBarry Smith for (l=0; l<nc; l++) { 106347c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 106447c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 106547c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1066aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 106747c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk); 106847c6ae99SBarry Smith } 106947c6ae99SBarry Smith } 107047c6ae99SBarry Smith } 107147c6ae99SBarry Smith } 107247c6ae99SBarry Smith rows[l] = l + nc*(slot); 107347c6ae99SBarry Smith } 107447c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 107547c6ae99SBarry Smith } 107647c6ae99SBarry Smith } 107747c6ae99SBarry Smith } 107847c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 107947c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 108047c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 108147c6ae99SBarry Smith } 108247c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 108347c6ae99SBarry Smith PetscFunctionReturn(0); 108447c6ae99SBarry Smith } 108547c6ae99SBarry Smith 108647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 108747c6ae99SBarry Smith 108847c6ae99SBarry Smith #undef __FUNCT__ 108994013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_1d_MPIAIJ" 109094013140SBarry Smith PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM da,Mat J) 109147c6ae99SBarry Smith { 109247c6ae99SBarry Smith PetscErrorCode ierr; 109347c6ae99SBarry Smith PetscInt xs,nx,i,i1,slot,gxs,gnx; 109447c6ae99SBarry Smith PetscInt m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l; 109547c6ae99SBarry Smith PetscInt istart,iend; 109647c6ae99SBarry Smith PetscScalar *values; 10971321219cSEthan Coon DMDABoundaryType bx; 109847c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 109947c6ae99SBarry Smith 110047c6ae99SBarry Smith PetscFunctionBegin; 110147c6ae99SBarry Smith /* 110247c6ae99SBarry Smith nc - number of components per grid point 110347c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 110447c6ae99SBarry Smith 110547c6ae99SBarry Smith */ 11061321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr); 110747c6ae99SBarry Smith col = 2*s + 1; 110847c6ae99SBarry Smith 1109aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr); 1110aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr); 111147c6ae99SBarry Smith 111247c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr); 111347c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr); 111447c6ae99SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 111547c6ae99SBarry Smith ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr); 111647c6ae99SBarry Smith 11171411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 11181411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 1119784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1120784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 112147c6ae99SBarry Smith 112247c6ae99SBarry Smith /* 112347c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 112447c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 112547c6ae99SBarry Smith PETSc ordering. 112647c6ae99SBarry Smith */ 1127fcfd50ebSBarry Smith if (!da->prealloc_only) { 112847c6ae99SBarry Smith ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 112947c6ae99SBarry Smith ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 113047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 113147c6ae99SBarry Smith istart = PetscMax(-s,gxs - i); 113247c6ae99SBarry Smith iend = PetscMin(s,gxs + gnx - i - 1); 113347c6ae99SBarry Smith slot = i - gxs; 113447c6ae99SBarry Smith 113547c6ae99SBarry Smith cnt = 0; 113647c6ae99SBarry Smith for (l=0; l<nc; l++) { 113747c6ae99SBarry Smith for (i1=istart; i1<iend+1; i1++) { 113847c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + i1); 113947c6ae99SBarry Smith } 114047c6ae99SBarry Smith rows[l] = l + nc*(slot); 114147c6ae99SBarry Smith } 114247c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 114347c6ae99SBarry Smith } 114447c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 114547c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 114647c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 114747c6ae99SBarry Smith } 114847c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 114947c6ae99SBarry Smith PetscFunctionReturn(0); 115047c6ae99SBarry Smith } 115147c6ae99SBarry Smith 115247c6ae99SBarry Smith #undef __FUNCT__ 115394013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIBAIJ" 115494013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM da,Mat J) 115547c6ae99SBarry Smith { 115647c6ae99SBarry Smith PetscErrorCode ierr; 115747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 115847c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz; 115947c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,ii,jj; 116047c6ae99SBarry Smith MPI_Comm comm; 116147c6ae99SBarry Smith PetscScalar *values; 11621321219cSEthan Coon DMDABoundaryType bx,by; 1163aa219208SBarry Smith DMDAStencilType st; 116447c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 116547c6ae99SBarry Smith 116647c6ae99SBarry Smith PetscFunctionBegin; 116747c6ae99SBarry Smith /* 116847c6ae99SBarry Smith nc - number of components per grid point 116947c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 117047c6ae99SBarry Smith */ 11711321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 117247c6ae99SBarry Smith col = 2*s + 1; 117347c6ae99SBarry Smith 1174aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 1175aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 117647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 117747c6ae99SBarry Smith 117847c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 117947c6ae99SBarry Smith 11801411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 11811411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 118247c6ae99SBarry Smith 118347c6ae99SBarry Smith /* determine the matrix preallocation information */ 118447c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr); 118547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 11861321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 11871321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 118847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 11891321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 11901321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 119147c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 119247c6ae99SBarry Smith 119347c6ae99SBarry Smith /* Find block columns in block row */ 119447c6ae99SBarry Smith cnt = 0; 119547c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 119647c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1197aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */ 119847c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 119947c6ae99SBarry Smith } 120047c6ae99SBarry Smith } 120147c6ae99SBarry Smith } 1202784ac674SJed Brown ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr); 120347c6ae99SBarry Smith } 120447c6ae99SBarry Smith } 120547c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 120647c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 120747c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 120847c6ae99SBarry Smith 1209784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1210784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 121147c6ae99SBarry Smith 121247c6ae99SBarry Smith /* 121347c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 121447c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 121547c6ae99SBarry Smith PETSc ordering. 121647c6ae99SBarry Smith */ 1217fcfd50ebSBarry Smith if (!da->prealloc_only) { 121847c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 121947c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 122047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 12211321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 12221321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 122347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 12241321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 12251321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 122647c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 122747c6ae99SBarry Smith cnt = 0; 122847c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 122947c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1230aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */ 123147c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 123247c6ae99SBarry Smith } 123347c6ae99SBarry Smith } 123447c6ae99SBarry Smith } 123547c6ae99SBarry Smith ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 123647c6ae99SBarry Smith } 123747c6ae99SBarry Smith } 123847c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 123947c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 124047c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 124147c6ae99SBarry Smith } 124247c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 124347c6ae99SBarry Smith PetscFunctionReturn(0); 124447c6ae99SBarry Smith } 124547c6ae99SBarry Smith 124647c6ae99SBarry Smith #undef __FUNCT__ 124794013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIBAIJ" 124894013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM da,Mat J) 124947c6ae99SBarry Smith { 125047c6ae99SBarry Smith PetscErrorCode ierr; 125147c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 125247c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz; 125347c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 125447c6ae99SBarry Smith MPI_Comm comm; 125547c6ae99SBarry Smith PetscScalar *values; 12561321219cSEthan Coon DMDABoundaryType bx,by,bz; 1257aa219208SBarry Smith DMDAStencilType st; 125847c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 125947c6ae99SBarry Smith 126047c6ae99SBarry Smith PetscFunctionBegin; 126147c6ae99SBarry Smith /* 126247c6ae99SBarry Smith nc - number of components per grid point 126347c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 126447c6ae99SBarry Smith 126547c6ae99SBarry Smith */ 12661321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 126747c6ae99SBarry Smith col = 2*s + 1; 126847c6ae99SBarry Smith 1269aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1270aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 127147c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 127247c6ae99SBarry Smith 127347c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr); 127447c6ae99SBarry Smith 12751411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 12761411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 127747c6ae99SBarry Smith 127847c6ae99SBarry Smith /* determine the matrix preallocation information */ 127947c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr); 128047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 12811321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 12821321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 128347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 12841321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 12851321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 128647c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 12871321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 12881321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 128947c6ae99SBarry Smith 129047c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 129147c6ae99SBarry Smith 129247c6ae99SBarry Smith /* Find block columns in block row */ 129347c6ae99SBarry Smith cnt = 0; 129447c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 129547c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 129647c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1297aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 129847c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 129947c6ae99SBarry Smith } 130047c6ae99SBarry Smith } 130147c6ae99SBarry Smith } 130247c6ae99SBarry Smith } 1303784ac674SJed Brown ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr); 130447c6ae99SBarry Smith } 130547c6ae99SBarry Smith } 130647c6ae99SBarry Smith } 130747c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 130847c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 130947c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 131047c6ae99SBarry Smith 1311784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1312784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 131347c6ae99SBarry Smith 131447c6ae99SBarry Smith /* 131547c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 131647c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 131747c6ae99SBarry Smith PETSc ordering. 131847c6ae99SBarry Smith */ 1319fcfd50ebSBarry Smith if (!da->prealloc_only) { 132047c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 132147c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 132247c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 13231321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 13241321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 132547c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 13261321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 13271321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 132847c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 13291321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 13301321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 133147c6ae99SBarry Smith 133247c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 133347c6ae99SBarry Smith 133447c6ae99SBarry Smith cnt = 0; 133547c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 133647c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 133747c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1338aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 133947c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 134047c6ae99SBarry Smith } 134147c6ae99SBarry Smith } 134247c6ae99SBarry Smith } 134347c6ae99SBarry Smith } 134447c6ae99SBarry Smith ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 134547c6ae99SBarry Smith } 134647c6ae99SBarry Smith } 134747c6ae99SBarry Smith } 134847c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 134947c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 135047c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 135147c6ae99SBarry Smith } 135247c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 135347c6ae99SBarry Smith PetscFunctionReturn(0); 135447c6ae99SBarry Smith } 135547c6ae99SBarry Smith 135647c6ae99SBarry Smith #undef __FUNCT__ 135747c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular" 135847c6ae99SBarry Smith /* 135947c6ae99SBarry Smith This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to 136047c6ae99SBarry Smith identify in the local ordering with periodic domain. 136147c6ae99SBarry Smith */ 136247c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[]) 136347c6ae99SBarry Smith { 136447c6ae99SBarry Smith PetscErrorCode ierr; 136547c6ae99SBarry Smith PetscInt i,n; 136647c6ae99SBarry Smith 136747c6ae99SBarry Smith PetscFunctionBegin; 136847c6ae99SBarry Smith ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr); 136947c6ae99SBarry Smith ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr); 137047c6ae99SBarry Smith for (i=0,n=0; i<*cnt; i++) { 137147c6ae99SBarry Smith if (col[i] >= *row) col[n++] = col[i]; 137247c6ae99SBarry Smith } 137347c6ae99SBarry Smith *cnt = n; 137447c6ae99SBarry Smith PetscFunctionReturn(0); 137547c6ae99SBarry Smith } 137647c6ae99SBarry Smith 137747c6ae99SBarry Smith #undef __FUNCT__ 137894013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPISBAIJ" 137994013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM da,Mat J) 138047c6ae99SBarry Smith { 138147c6ae99SBarry Smith PetscErrorCode ierr; 138247c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 138347c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz; 138447c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,ii,jj; 138547c6ae99SBarry Smith MPI_Comm comm; 138647c6ae99SBarry Smith PetscScalar *values; 13871321219cSEthan Coon DMDABoundaryType bx,by; 1388aa219208SBarry Smith DMDAStencilType st; 138947c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 139047c6ae99SBarry Smith 139147c6ae99SBarry Smith PetscFunctionBegin; 139247c6ae99SBarry Smith /* 139347c6ae99SBarry Smith nc - number of components per grid point 139447c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 139547c6ae99SBarry Smith */ 13961321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 139747c6ae99SBarry Smith col = 2*s + 1; 139847c6ae99SBarry Smith 1399aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 1400aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 140147c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 140247c6ae99SBarry Smith 140347c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 140447c6ae99SBarry Smith 14051411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 14061411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 140747c6ae99SBarry Smith 140847c6ae99SBarry Smith /* determine the matrix preallocation information */ 140947c6ae99SBarry Smith ierr = MatPreallocateSymmetricInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr); 141047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 14111321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 14121321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 141347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 14141321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 14151321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 141647c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 141747c6ae99SBarry Smith 141847c6ae99SBarry Smith /* Find block columns in block row */ 141947c6ae99SBarry Smith cnt = 0; 142047c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 142147c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1422aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { 142347c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 142447c6ae99SBarry Smith } 142547c6ae99SBarry Smith } 142647c6ae99SBarry Smith } 142747c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 142847c6ae99SBarry Smith ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr); 142947c6ae99SBarry Smith } 143047c6ae99SBarry Smith } 143147c6ae99SBarry Smith ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 143247c6ae99SBarry Smith ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 143347c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 143447c6ae99SBarry Smith 1435784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1436784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 143747c6ae99SBarry Smith 143847c6ae99SBarry Smith /* 143947c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 144047c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 144147c6ae99SBarry Smith PETSc ordering. 144247c6ae99SBarry Smith */ 1443fcfd50ebSBarry Smith if (!da->prealloc_only) { 144447c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 144547c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 144647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 14471321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 14481321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 144947c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 14501321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 14511321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 145247c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 145347c6ae99SBarry Smith 145447c6ae99SBarry Smith /* Find block columns in block row */ 145547c6ae99SBarry Smith cnt = 0; 145647c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 145747c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1458aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { 145947c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 146047c6ae99SBarry Smith } 146147c6ae99SBarry Smith } 146247c6ae99SBarry Smith } 146347c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 146447c6ae99SBarry Smith ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 146547c6ae99SBarry Smith } 146647c6ae99SBarry Smith } 146747c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 146847c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 146947c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147047c6ae99SBarry Smith } 147147c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 147247c6ae99SBarry Smith PetscFunctionReturn(0); 147347c6ae99SBarry Smith } 147447c6ae99SBarry Smith 147547c6ae99SBarry Smith #undef __FUNCT__ 147694013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPISBAIJ" 147794013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPISBAIJ(DM da,Mat J) 147847c6ae99SBarry Smith { 147947c6ae99SBarry Smith PetscErrorCode ierr; 148047c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 148147c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz; 148247c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 148347c6ae99SBarry Smith MPI_Comm comm; 148447c6ae99SBarry Smith PetscScalar *values; 14851321219cSEthan Coon DMDABoundaryType bx,by,bz; 1486aa219208SBarry Smith DMDAStencilType st; 148747c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 148847c6ae99SBarry Smith 148947c6ae99SBarry Smith PetscFunctionBegin; 149047c6ae99SBarry Smith /* 149147c6ae99SBarry Smith nc - number of components per grid point 149247c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 149347c6ae99SBarry Smith */ 14941321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 149547c6ae99SBarry Smith col = 2*s + 1; 149647c6ae99SBarry Smith 1497aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1498aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 149947c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 150047c6ae99SBarry Smith 150147c6ae99SBarry Smith /* create the matrix */ 150247c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr); 150347c6ae99SBarry Smith 15041411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 15051411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 150647c6ae99SBarry Smith 150747c6ae99SBarry Smith /* determine the matrix preallocation information */ 150847c6ae99SBarry Smith ierr = MatPreallocateSymmetricInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr); 150947c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 15101321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 15111321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 151247c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 15131321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 15141321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 151547c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 15161321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 15171321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 151847c6ae99SBarry Smith 151947c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 152047c6ae99SBarry Smith 152147c6ae99SBarry Smith /* Find block columns in block row */ 152247c6ae99SBarry Smith cnt = 0; 152347c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 152447c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 152547c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1526aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) { 152747c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 152847c6ae99SBarry Smith } 152947c6ae99SBarry Smith } 153047c6ae99SBarry Smith } 153147c6ae99SBarry Smith } 153247c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 153347c6ae99SBarry Smith ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr); 153447c6ae99SBarry Smith } 153547c6ae99SBarry Smith } 153647c6ae99SBarry Smith } 153747c6ae99SBarry Smith ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 153847c6ae99SBarry Smith ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 153947c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 154047c6ae99SBarry Smith 1541784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1542784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 154347c6ae99SBarry Smith 154447c6ae99SBarry Smith /* 154547c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 154647c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 154747c6ae99SBarry Smith PETSc ordering. 154847c6ae99SBarry Smith */ 1549fcfd50ebSBarry Smith if (!da->prealloc_only) { 155047c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 155147c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 155247c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 15531321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 15541321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 155547c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 15561321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 15571321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 155847c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 15591321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 15601321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 156147c6ae99SBarry Smith 156247c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 156347c6ae99SBarry Smith 156447c6ae99SBarry Smith cnt = 0; 156547c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 156647c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 156747c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1568aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) { 156947c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 157047c6ae99SBarry Smith } 157147c6ae99SBarry Smith } 157247c6ae99SBarry Smith } 157347c6ae99SBarry Smith } 157447c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 157547c6ae99SBarry Smith ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 157647c6ae99SBarry Smith } 157747c6ae99SBarry Smith } 157847c6ae99SBarry Smith } 157947c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 158047c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 158147c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 158247c6ae99SBarry Smith } 158347c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 158447c6ae99SBarry Smith PetscFunctionReturn(0); 158547c6ae99SBarry Smith } 158647c6ae99SBarry Smith 158747c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 158847c6ae99SBarry Smith 158947c6ae99SBarry Smith #undef __FUNCT__ 159094013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIAIJ_Fill" 159194013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J) 159247c6ae99SBarry Smith { 159347c6ae99SBarry Smith PetscErrorCode ierr; 159447c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 159547c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz; 159647c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 159747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 159847c6ae99SBarry Smith PetscInt ifill_col,*dfill = dd->dfill,*ofill = dd->ofill; 159947c6ae99SBarry Smith MPI_Comm comm; 160047c6ae99SBarry Smith PetscScalar *values; 16011321219cSEthan Coon DMDABoundaryType bx,by,bz; 160247c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 1603aa219208SBarry Smith DMDAStencilType st; 160447c6ae99SBarry Smith 160547c6ae99SBarry Smith PetscFunctionBegin; 160647c6ae99SBarry Smith /* 160747c6ae99SBarry Smith nc - number of components per grid point 160847c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 160947c6ae99SBarry Smith 161047c6ae99SBarry Smith */ 16111321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 161247c6ae99SBarry Smith col = 2*s + 1; 16131321219cSEthan Coon if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){ 161447c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\ 161547c6ae99SBarry Smith by 2*stencil_width + 1\n"); 161647c6ae99SBarry Smith } 16171321219cSEthan Coon if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){ 161847c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\ 161947c6ae99SBarry Smith by 2*stencil_width + 1\n"); 162047c6ae99SBarry Smith } 16211321219cSEthan Coon if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)){ 162247c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\ 162347c6ae99SBarry Smith by 2*stencil_width + 1\n"); 162447c6ae99SBarry Smith } 162547c6ae99SBarry Smith 1626aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1627aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 162847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 162947c6ae99SBarry Smith 163047c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 16311411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 16321411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 163347c6ae99SBarry Smith 163447c6ae99SBarry Smith /* determine the matrix preallocation information */ 163547c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr); 163647c6ae99SBarry Smith 163747c6ae99SBarry Smith 163847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 16391321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 16401321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 164147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 16421321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 16431321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 164447c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 16451321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 16461321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 164747c6ae99SBarry Smith 164847c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 164947c6ae99SBarry Smith 165047c6ae99SBarry Smith for (l=0; l<nc; l++) { 165147c6ae99SBarry Smith cnt = 0; 165247c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 165347c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 165447c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 165547c6ae99SBarry Smith if (ii || jj || kk) { 1656aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 165747c6ae99SBarry Smith for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) 165847c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 165947c6ae99SBarry Smith } 166047c6ae99SBarry Smith } else { 166147c6ae99SBarry Smith if (dfill) { 166247c6ae99SBarry Smith for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) 166347c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 166447c6ae99SBarry Smith } else { 166547c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 166647c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk); 166747c6ae99SBarry Smith } 166847c6ae99SBarry Smith } 166947c6ae99SBarry Smith } 167047c6ae99SBarry Smith } 167147c6ae99SBarry Smith } 167247c6ae99SBarry Smith row = l + nc*(slot); 1673784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 167447c6ae99SBarry Smith } 167547c6ae99SBarry Smith } 167647c6ae99SBarry Smith } 167747c6ae99SBarry Smith } 167847c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 167947c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 168047c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1681784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1682784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 168347c6ae99SBarry Smith 168447c6ae99SBarry Smith /* 168547c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 168647c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 168747c6ae99SBarry Smith PETSc ordering. 168847c6ae99SBarry Smith */ 1689fcfd50ebSBarry Smith if (!da->prealloc_only) { 169047c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 169147c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 169247c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 16931321219cSEthan Coon istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 16941321219cSEthan Coon iend = (bx == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 169547c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 16961321219cSEthan Coon jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 16971321219cSEthan Coon jend = (by == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 169847c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 16991321219cSEthan Coon kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 17001321219cSEthan Coon kend = (bz == DMDA_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 170147c6ae99SBarry Smith 170247c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 170347c6ae99SBarry Smith 170447c6ae99SBarry Smith for (l=0; l<nc; l++) { 170547c6ae99SBarry Smith cnt = 0; 170647c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 170747c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 170847c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 170947c6ae99SBarry Smith if (ii || jj || kk) { 1710aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 171147c6ae99SBarry Smith for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) 171247c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 171347c6ae99SBarry Smith } 171447c6ae99SBarry Smith } else { 171547c6ae99SBarry Smith if (dfill) { 171647c6ae99SBarry Smith for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) 171747c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 171847c6ae99SBarry Smith } else { 171947c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 172047c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk); 172147c6ae99SBarry Smith } 172247c6ae99SBarry Smith } 172347c6ae99SBarry Smith } 172447c6ae99SBarry Smith } 172547c6ae99SBarry Smith } 172647c6ae99SBarry Smith row = l + nc*(slot); 172747c6ae99SBarry Smith ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 172847c6ae99SBarry Smith } 172947c6ae99SBarry Smith } 173047c6ae99SBarry Smith } 173147c6ae99SBarry Smith } 173247c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 173347c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 173447c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 173547c6ae99SBarry Smith } 173647c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 173747c6ae99SBarry Smith PetscFunctionReturn(0); 173847c6ae99SBarry Smith } 1739