147c6ae99SBarry Smith #define PETSCDM_DLL 247c6ae99SBarry Smith 3e1589f56SBarry Smith #include "private/daimpl.h" /*I "petscdm.h" I*/ 447c6ae99SBarry Smith #include "petscmat.h" /*I "petscmat.h" I*/ 547c6ae99SBarry Smith #include "private/matimpl.h" 647c6ae99SBarry Smith 7*09573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring *); 8*09573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring *); 9*09573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring *); 10*09573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_3d_MPIAIJ(DM,ISColoringType,ISColoring *); 1147c6ae99SBarry Smith 1247c6ae99SBarry Smith /* 1347c6ae99SBarry Smith For ghost i that may be negative or greater than the upper bound this 1447c6ae99SBarry Smith maps it into the 0:m-1 range using periodicity 1547c6ae99SBarry Smith */ 1647c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i:((i >= m) ? i-m:i)) 1747c6ae99SBarry Smith 1847c6ae99SBarry Smith #undef __FUNCT__ 19aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills_Private" 20aa219208SBarry Smith static PetscErrorCode DMDASetBlockFills_Private(PetscInt *dfill,PetscInt w,PetscInt **rfill) 2147c6ae99SBarry Smith { 2247c6ae99SBarry Smith PetscErrorCode ierr; 2347c6ae99SBarry Smith PetscInt i,j,nz,*fill; 2447c6ae99SBarry Smith 2547c6ae99SBarry Smith PetscFunctionBegin; 2647c6ae99SBarry Smith if (!dfill) PetscFunctionReturn(0); 2747c6ae99SBarry Smith 2847c6ae99SBarry Smith /* count number nonzeros */ 2947c6ae99SBarry Smith nz = 0; 3047c6ae99SBarry Smith for (i=0; i<w; i++) { 3147c6ae99SBarry Smith for (j=0; j<w; j++) { 3247c6ae99SBarry Smith if (dfill[w*i+j]) nz++; 3347c6ae99SBarry Smith } 3447c6ae99SBarry Smith } 3547c6ae99SBarry Smith ierr = PetscMalloc((nz + w + 1)*sizeof(PetscInt),&fill);CHKERRQ(ierr); 3647c6ae99SBarry Smith /* construct modified CSR storage of nonzero structure */ 3747c6ae99SBarry Smith nz = w + 1; 3847c6ae99SBarry Smith for (i=0; i<w; i++) { 3947c6ae99SBarry Smith fill[i] = nz; 4047c6ae99SBarry Smith for (j=0; j<w; j++) { 4147c6ae99SBarry Smith if (dfill[w*i+j]) { 4247c6ae99SBarry Smith fill[nz] = j; 4347c6ae99SBarry Smith nz++; 4447c6ae99SBarry Smith } 4547c6ae99SBarry Smith } 4647c6ae99SBarry Smith } 4747c6ae99SBarry Smith fill[w] = nz; 4847c6ae99SBarry Smith 4947c6ae99SBarry Smith *rfill = fill; 5047c6ae99SBarry Smith PetscFunctionReturn(0); 5147c6ae99SBarry Smith } 5247c6ae99SBarry Smith 5347c6ae99SBarry Smith #undef __FUNCT__ 54aa219208SBarry Smith #define __FUNCT__ "DMDASetMatPreallocateOnly" 5547c6ae99SBarry Smith /*@ 56aa219208SBarry Smith DMDASetMatPreallocateOnly - When DMGetMatrix() is called the matrix will be properly 5747c6ae99SBarry Smith preallocated but the nonzero structure and zero values will not be set. 5847c6ae99SBarry Smith 59aa219208SBarry Smith Logically Collective on DMDA 6047c6ae99SBarry Smith 6147c6ae99SBarry Smith Input Parameter: 6247c6ae99SBarry Smith + da - the distributed array 6347c6ae99SBarry Smith - only - PETSC_TRUE if only want preallocation 6447c6ae99SBarry Smith 6547c6ae99SBarry Smith 6647c6ae99SBarry Smith Level: developer 6747c6ae99SBarry Smith 68aa219208SBarry Smith .seealso DMGetMatrix(), DMDASetGetMatrix(), DMDASetBlockSize(), DMDASetBlockFills() 6947c6ae99SBarry Smith 7047c6ae99SBarry Smith @*/ 71aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetMatPreallocateOnly(DM da,PetscBool only) 7247c6ae99SBarry Smith { 7347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 7447c6ae99SBarry Smith 7547c6ae99SBarry Smith PetscFunctionBegin; 7647c6ae99SBarry Smith dd->prealloc_only = only; 7747c6ae99SBarry Smith PetscFunctionReturn(0); 7847c6ae99SBarry Smith } 7947c6ae99SBarry Smith 8047c6ae99SBarry Smith #undef __FUNCT__ 81aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills" 8247c6ae99SBarry Smith /*@ 83aa219208SBarry Smith DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem 8494013140SBarry Smith of the matrix returned by DMGetMatrix(). 8547c6ae99SBarry Smith 86aa219208SBarry Smith Logically Collective on DMDA 8747c6ae99SBarry Smith 8847c6ae99SBarry Smith Input Parameter: 8947c6ae99SBarry Smith + da - the distributed array 9047c6ae99SBarry Smith . dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block) 9147c6ae99SBarry Smith - ofill - the fill pattern in the off-diagonal blocks 9247c6ae99SBarry Smith 9347c6ae99SBarry Smith 9447c6ae99SBarry Smith Level: developer 9547c6ae99SBarry Smith 9647c6ae99SBarry Smith Notes: This only makes sense when you are doing multicomponent problems but using the 9747c6ae99SBarry Smith MPIAIJ matrix format 9847c6ae99SBarry Smith 9947c6ae99SBarry Smith The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries 10047c6ae99SBarry Smith representing coupling and 0 entries for missing coupling. For example 10147c6ae99SBarry Smith $ dfill[9] = {1, 0, 0, 10247c6ae99SBarry Smith $ 1, 1, 0, 10347c6ae99SBarry Smith $ 0, 1, 1} 10447c6ae99SBarry Smith means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with 10547c6ae99SBarry Smith itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the 10647c6ae99SBarry Smith diagonal block). 10747c6ae99SBarry Smith 108aa219208SBarry Smith DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then 10947c6ae99SBarry Smith can be represented in the dfill, ofill format 11047c6ae99SBarry Smith 11147c6ae99SBarry Smith Contributed by Glenn Hammond 11247c6ae99SBarry Smith 113aa219208SBarry Smith .seealso DMGetMatrix(), DMDASetGetMatrix(), DMDASetMatPreallocateOnly() 11447c6ae99SBarry Smith 11547c6ae99SBarry Smith @*/ 116aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetBlockFills(DM da,PetscInt *dfill,PetscInt *ofill) 11747c6ae99SBarry Smith { 11847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 11947c6ae99SBarry Smith PetscErrorCode ierr; 12047c6ae99SBarry Smith 12147c6ae99SBarry Smith PetscFunctionBegin; 122aa219208SBarry Smith ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr); 123aa219208SBarry Smith ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr); 12447c6ae99SBarry Smith PetscFunctionReturn(0); 12547c6ae99SBarry Smith } 12647c6ae99SBarry Smith 12747c6ae99SBarry Smith 12847c6ae99SBarry Smith #undef __FUNCT__ 12994013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA" 13094013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring_DA(DM da,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 13147c6ae99SBarry Smith { 13247c6ae99SBarry Smith PetscErrorCode ierr; 13347c6ae99SBarry Smith PetscInt dim,m,n,p,nc; 134aa219208SBarry Smith DMDAPeriodicType wrap; 13547c6ae99SBarry Smith MPI_Comm comm; 13647c6ae99SBarry Smith PetscMPIInt size; 13747c6ae99SBarry Smith PetscBool isBAIJ; 13847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 13947c6ae99SBarry Smith 14047c6ae99SBarry Smith PetscFunctionBegin; 14147c6ae99SBarry Smith /* 14247c6ae99SBarry Smith m 14347c6ae99SBarry Smith ------------------------------------------------------ 14447c6ae99SBarry Smith | | 14547c6ae99SBarry Smith | | 14647c6ae99SBarry Smith | ---------------------- | 14747c6ae99SBarry Smith | | | | 14847c6ae99SBarry Smith n | yn | | | 14947c6ae99SBarry Smith | | | | 15047c6ae99SBarry Smith | .--------------------- | 15147c6ae99SBarry Smith | (xs,ys) xn | 15247c6ae99SBarry Smith | . | 15347c6ae99SBarry Smith | (gxs,gys) | 15447c6ae99SBarry Smith | | 15547c6ae99SBarry Smith ----------------------------------------------------- 15647c6ae99SBarry Smith */ 15747c6ae99SBarry Smith 15847c6ae99SBarry Smith /* 15947c6ae99SBarry Smith nc - number of components per grid point 16047c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 16147c6ae99SBarry Smith 16247c6ae99SBarry Smith */ 163aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&wrap,0);CHKERRQ(ierr); 16447c6ae99SBarry Smith 16547c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 16647c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 16747c6ae99SBarry Smith if (ctype == IS_COLORING_GHOSTED){ 16847c6ae99SBarry Smith if (size == 1) { 16947c6ae99SBarry Smith ctype = IS_COLORING_GLOBAL; 17047c6ae99SBarry Smith } else if (dim > 1){ 171aa219208SBarry Smith if ((m==1 && DMDAXPeriodic(wrap)) || (n==1 && DMDAYPeriodic(wrap)) || (p==1 && DMDAZPeriodic(wrap))){ 17247c6ae99SBarry 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"); 17347c6ae99SBarry Smith } 17447c6ae99SBarry Smith } 17547c6ae99SBarry Smith } 17647c6ae99SBarry Smith 177aa219208SBarry Smith /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ 17847c6ae99SBarry Smith matrices is for the blocks, not the individual matrix elements */ 17947c6ae99SBarry Smith ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr); 18047c6ae99SBarry Smith if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);} 18147c6ae99SBarry Smith if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);} 18247c6ae99SBarry Smith if (isBAIJ) { 18347c6ae99SBarry Smith dd->w = 1; 18447c6ae99SBarry Smith dd->xs = dd->xs/nc; 18547c6ae99SBarry Smith dd->xe = dd->xe/nc; 18647c6ae99SBarry Smith dd->Xs = dd->Xs/nc; 18747c6ae99SBarry Smith dd->Xe = dd->Xe/nc; 18847c6ae99SBarry Smith } 18947c6ae99SBarry Smith 19047c6ae99SBarry Smith /* 191aa219208SBarry Smith We do not provide a getcoloring function in the DMDA operations because 192aa219208SBarry Smith the basic DMDA does not know about matrices. We think of DMDA as being more 19347c6ae99SBarry Smith more low-level then matrices. 19447c6ae99SBarry Smith */ 19547c6ae99SBarry Smith if (dim == 1) { 19694013140SBarry Smith ierr = DMGetColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 19747c6ae99SBarry Smith } else if (dim == 2) { 19894013140SBarry Smith ierr = DMGetColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 19947c6ae99SBarry Smith } else if (dim == 3) { 20094013140SBarry Smith ierr = DMGetColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 20147c6ae99SBarry Smith } else { 20247c6ae99SBarry Smith SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim); 20347c6ae99SBarry Smith } 20447c6ae99SBarry Smith if (isBAIJ) { 20547c6ae99SBarry Smith dd->w = nc; 20647c6ae99SBarry Smith dd->xs = dd->xs*nc; 20747c6ae99SBarry Smith dd->xe = dd->xe*nc; 20847c6ae99SBarry Smith dd->Xs = dd->Xs*nc; 20947c6ae99SBarry Smith dd->Xe = dd->Xe*nc; 21047c6ae99SBarry Smith } 21147c6ae99SBarry Smith PetscFunctionReturn(0); 21247c6ae99SBarry Smith } 21347c6ae99SBarry Smith 21447c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 21547c6ae99SBarry Smith 21647c6ae99SBarry Smith #undef __FUNCT__ 21794013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_2d_MPIAIJ" 21894013140SBarry Smith PetscErrorCode DMGetColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 21947c6ae99SBarry Smith { 22047c6ae99SBarry Smith PetscErrorCode ierr; 22147c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col; 22247c6ae99SBarry Smith PetscInt ncolors; 22347c6ae99SBarry Smith MPI_Comm comm; 224aa219208SBarry Smith DMDAPeriodicType wrap; 225aa219208SBarry Smith DMDAStencilType st; 22647c6ae99SBarry Smith ISColoringValue *colors; 22747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 22847c6ae99SBarry Smith 22947c6ae99SBarry Smith PetscFunctionBegin; 23047c6ae99SBarry Smith /* 23147c6ae99SBarry Smith nc - number of components per grid point 23247c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 23347c6ae99SBarry Smith 23447c6ae99SBarry Smith */ 235aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&wrap,&st);CHKERRQ(ierr); 23647c6ae99SBarry Smith col = 2*s + 1; 237aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 238aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 23947c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 24047c6ae99SBarry Smith 24147c6ae99SBarry Smith /* special case as taught to us by Paul Hovland */ 242aa219208SBarry Smith if (st == DMDA_STENCIL_STAR && s == 1) { 24394013140SBarry Smith ierr = DMGetColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 24447c6ae99SBarry Smith } else { 24547c6ae99SBarry Smith 246aa219208SBarry Smith if (DMDAXPeriodic(wrap) && (m % col)){ 24747c6ae99SBarry Smith SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\ 24847c6ae99SBarry Smith by 2*stencil_width + 1 (%d)\n", m, col); 24947c6ae99SBarry Smith } 250aa219208SBarry Smith if (DMDAYPeriodic(wrap) && (n % col)){ 25147c6ae99SBarry Smith SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\ 25247c6ae99SBarry Smith by 2*stencil_width + 1 (%d)\n", n, col); 25347c6ae99SBarry Smith } 25447c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 25547c6ae99SBarry Smith if (!dd->localcoloring) { 25647c6ae99SBarry Smith ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 25747c6ae99SBarry Smith ii = 0; 25847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 25947c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 26047c6ae99SBarry Smith for (k=0; k<nc; k++) { 26147c6ae99SBarry Smith colors[ii++] = k + nc*((i % col) + col*(j % col)); 26247c6ae99SBarry Smith } 26347c6ae99SBarry Smith } 26447c6ae99SBarry Smith } 26547c6ae99SBarry Smith ncolors = nc + nc*(col-1 + col*(col-1)); 26647c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr); 26747c6ae99SBarry Smith } 26847c6ae99SBarry Smith *coloring = dd->localcoloring; 26947c6ae99SBarry Smith } else if (ctype == IS_COLORING_GHOSTED) { 27047c6ae99SBarry Smith if (!dd->ghostedcoloring) { 27147c6ae99SBarry Smith ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 27247c6ae99SBarry Smith ii = 0; 27347c6ae99SBarry Smith for (j=gys; j<gys+gny; j++) { 27447c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 27547c6ae99SBarry Smith for (k=0; k<nc; k++) { 27647c6ae99SBarry Smith /* the complicated stuff is to handle periodic boundaries */ 27747c6ae99SBarry Smith colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col)); 27847c6ae99SBarry Smith } 27947c6ae99SBarry Smith } 28047c6ae99SBarry Smith } 28147c6ae99SBarry Smith ncolors = nc + nc*(col - 1 + col*(col-1)); 28247c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr); 28347c6ae99SBarry Smith /* PetscIntView(ncolors,(PetscInt *)colors,0); */ 28447c6ae99SBarry Smith 28547c6ae99SBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr); 28647c6ae99SBarry Smith } 28747c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 28847c6ae99SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 28947c6ae99SBarry Smith } 29047c6ae99SBarry Smith ierr = ISColoringReference(*coloring);CHKERRQ(ierr); 29147c6ae99SBarry Smith PetscFunctionReturn(0); 29247c6ae99SBarry Smith } 29347c6ae99SBarry Smith 29447c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 29547c6ae99SBarry Smith 29647c6ae99SBarry Smith #undef __FUNCT__ 29794013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_3d_MPIAIJ" 29894013140SBarry Smith PetscErrorCode DMGetColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 29947c6ae99SBarry Smith { 30047c6ae99SBarry Smith PetscErrorCode ierr; 30147c6ae99SBarry 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; 30247c6ae99SBarry Smith PetscInt ncolors; 30347c6ae99SBarry Smith MPI_Comm comm; 304aa219208SBarry Smith DMDAPeriodicType wrap; 305aa219208SBarry Smith DMDAStencilType st; 30647c6ae99SBarry Smith ISColoringValue *colors; 30747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 30847c6ae99SBarry Smith 30947c6ae99SBarry Smith PetscFunctionBegin; 31047c6ae99SBarry Smith /* 31147c6ae99SBarry Smith nc - number of components per grid point 31247c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 31347c6ae99SBarry Smith 31447c6ae99SBarry Smith */ 315aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&wrap,&st);CHKERRQ(ierr); 31647c6ae99SBarry Smith col = 2*s + 1; 317aa219208SBarry Smith if (DMDAXPeriodic(wrap) && (m % col)){ 31847c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\ 31947c6ae99SBarry Smith by 2*stencil_width + 1\n"); 32047c6ae99SBarry Smith } 321aa219208SBarry Smith if (DMDAYPeriodic(wrap) && (n % col)){ 32247c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\ 32347c6ae99SBarry Smith by 2*stencil_width + 1\n"); 32447c6ae99SBarry Smith } 325aa219208SBarry Smith if (DMDAZPeriodic(wrap) && (p % col)){ 32647c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\ 32747c6ae99SBarry Smith by 2*stencil_width + 1\n"); 32847c6ae99SBarry Smith } 32947c6ae99SBarry Smith 330aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 331aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 33247c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 33347c6ae99SBarry Smith 33447c6ae99SBarry Smith /* create the coloring */ 33547c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 33647c6ae99SBarry Smith if (!dd->localcoloring) { 33747c6ae99SBarry Smith ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 33847c6ae99SBarry Smith ii = 0; 33947c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 34047c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 34147c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 34247c6ae99SBarry Smith for (l=0; l<nc; l++) { 34347c6ae99SBarry Smith colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col)); 34447c6ae99SBarry Smith } 34547c6ae99SBarry Smith } 34647c6ae99SBarry Smith } 34747c6ae99SBarry Smith } 34847c6ae99SBarry Smith ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1)); 34947c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr); 35047c6ae99SBarry Smith } 35147c6ae99SBarry Smith *coloring = dd->localcoloring; 35247c6ae99SBarry Smith } else if (ctype == IS_COLORING_GHOSTED) { 35347c6ae99SBarry Smith if (!dd->ghostedcoloring) { 35447c6ae99SBarry Smith ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 35547c6ae99SBarry Smith ii = 0; 35647c6ae99SBarry Smith for (k=gzs; k<gzs+gnz; k++) { 35747c6ae99SBarry Smith for (j=gys; j<gys+gny; j++) { 35847c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 35947c6ae99SBarry Smith for (l=0; l<nc; l++) { 36047c6ae99SBarry Smith /* the complicated stuff is to handle periodic boundaries */ 36147c6ae99SBarry Smith colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col)); 36247c6ae99SBarry Smith } 36347c6ae99SBarry Smith } 36447c6ae99SBarry Smith } 36547c6ae99SBarry Smith } 36647c6ae99SBarry Smith ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1)); 36747c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr); 36847c6ae99SBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr); 36947c6ae99SBarry Smith } 37047c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 37147c6ae99SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 37247c6ae99SBarry Smith ierr = ISColoringReference(*coloring);CHKERRQ(ierr); 37347c6ae99SBarry Smith PetscFunctionReturn(0); 37447c6ae99SBarry Smith } 37547c6ae99SBarry Smith 37647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 37747c6ae99SBarry Smith 37847c6ae99SBarry Smith #undef __FUNCT__ 37994013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_1d_MPIAIJ" 38094013140SBarry Smith PetscErrorCode DMGetColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 38147c6ae99SBarry Smith { 38247c6ae99SBarry Smith PetscErrorCode ierr; 38347c6ae99SBarry Smith PetscInt xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col; 38447c6ae99SBarry Smith PetscInt ncolors; 38547c6ae99SBarry Smith MPI_Comm comm; 386aa219208SBarry Smith DMDAPeriodicType wrap; 38747c6ae99SBarry Smith ISColoringValue *colors; 38847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 38947c6ae99SBarry Smith 39047c6ae99SBarry Smith PetscFunctionBegin; 39147c6ae99SBarry Smith /* 39247c6ae99SBarry Smith nc - number of components per grid point 39347c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 39447c6ae99SBarry Smith 39547c6ae99SBarry Smith */ 396aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&wrap,0);CHKERRQ(ierr); 39747c6ae99SBarry Smith col = 2*s + 1; 39847c6ae99SBarry Smith 399aa219208SBarry Smith if (DMDAXPeriodic(wrap) && (m % col)) { 40047c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points is divisible\n\ 40147c6ae99SBarry Smith by 2*stencil_width + 1\n"); 40247c6ae99SBarry Smith } 40347c6ae99SBarry Smith 404aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr); 405aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr); 40647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 40747c6ae99SBarry Smith 40847c6ae99SBarry Smith /* create the coloring */ 40947c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 41047c6ae99SBarry Smith if (!dd->localcoloring) { 41147c6ae99SBarry Smith ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 41247c6ae99SBarry Smith i1 = 0; 41347c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 41447c6ae99SBarry Smith for (l=0; l<nc; l++) { 41547c6ae99SBarry Smith colors[i1++] = l + nc*(i % col); 41647c6ae99SBarry Smith } 41747c6ae99SBarry Smith } 41847c6ae99SBarry Smith ncolors = nc + nc*(col-1); 41947c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr); 42047c6ae99SBarry Smith } 42147c6ae99SBarry Smith *coloring = dd->localcoloring; 42247c6ae99SBarry Smith } else if (ctype == IS_COLORING_GHOSTED) { 42347c6ae99SBarry Smith if (!dd->ghostedcoloring) { 42447c6ae99SBarry Smith ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 42547c6ae99SBarry Smith i1 = 0; 42647c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 42747c6ae99SBarry Smith for (l=0; l<nc; l++) { 42847c6ae99SBarry Smith /* the complicated stuff is to handle periodic boundaries */ 42947c6ae99SBarry Smith colors[i1++] = l + nc*(SetInRange(i,m) % col); 43047c6ae99SBarry Smith } 43147c6ae99SBarry Smith } 43247c6ae99SBarry Smith ncolors = nc + nc*(col-1); 43347c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr); 43447c6ae99SBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr); 43547c6ae99SBarry Smith } 43647c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 43747c6ae99SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 43847c6ae99SBarry Smith ierr = ISColoringReference(*coloring);CHKERRQ(ierr); 43947c6ae99SBarry Smith PetscFunctionReturn(0); 44047c6ae99SBarry Smith } 44147c6ae99SBarry Smith 44247c6ae99SBarry Smith #undef __FUNCT__ 44394013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_2d_5pt_MPIAIJ" 44494013140SBarry Smith PetscErrorCode DMGetColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 44547c6ae99SBarry Smith { 44647c6ae99SBarry Smith PetscErrorCode ierr; 44747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc; 44847c6ae99SBarry Smith PetscInt ncolors; 44947c6ae99SBarry Smith MPI_Comm comm; 450aa219208SBarry Smith DMDAPeriodicType wrap; 45147c6ae99SBarry Smith ISColoringValue *colors; 45247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 45347c6ae99SBarry Smith 45447c6ae99SBarry Smith PetscFunctionBegin; 45547c6ae99SBarry Smith /* 45647c6ae99SBarry Smith nc - number of components per grid point 45747c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 45847c6ae99SBarry Smith 45947c6ae99SBarry Smith */ 460aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,0);CHKERRQ(ierr); 461aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 462aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 46347c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 46447c6ae99SBarry Smith 465aa219208SBarry Smith if (DMDAXPeriodic(wrap) && (m % 5)){ 46647c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\ 46747c6ae99SBarry Smith by 5\n"); 46847c6ae99SBarry Smith } 469aa219208SBarry Smith if (DMDAYPeriodic(wrap) && (n % 5)){ 47047c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\ 47147c6ae99SBarry Smith by 5\n"); 47247c6ae99SBarry Smith } 47347c6ae99SBarry Smith 47447c6ae99SBarry Smith /* create the coloring */ 47547c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 47647c6ae99SBarry Smith if (!dd->localcoloring) { 47747c6ae99SBarry Smith ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 47847c6ae99SBarry Smith ii = 0; 47947c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 48047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 48147c6ae99SBarry Smith for (k=0; k<nc; k++) { 48247c6ae99SBarry Smith colors[ii++] = k + nc*((3*j+i) % 5); 48347c6ae99SBarry Smith } 48447c6ae99SBarry Smith } 48547c6ae99SBarry Smith } 48647c6ae99SBarry Smith ncolors = 5*nc; 48747c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr); 48847c6ae99SBarry Smith } 48947c6ae99SBarry Smith *coloring = dd->localcoloring; 49047c6ae99SBarry Smith } else if (ctype == IS_COLORING_GHOSTED) { 49147c6ae99SBarry Smith if (!dd->ghostedcoloring) { 49247c6ae99SBarry Smith ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); 49347c6ae99SBarry Smith ii = 0; 49447c6ae99SBarry Smith for (j=gys; j<gys+gny; j++) { 49547c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 49647c6ae99SBarry Smith for (k=0; k<nc; k++) { 49747c6ae99SBarry Smith colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5); 49847c6ae99SBarry Smith } 49947c6ae99SBarry Smith } 50047c6ae99SBarry Smith } 50147c6ae99SBarry Smith ncolors = 5*nc; 50247c6ae99SBarry Smith ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr); 50347c6ae99SBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr); 50447c6ae99SBarry Smith } 50547c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 50647c6ae99SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 50747c6ae99SBarry Smith PetscFunctionReturn(0); 50847c6ae99SBarry Smith } 50947c6ae99SBarry Smith 51047c6ae99SBarry Smith /* =========================================================================== */ 511*09573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM,Mat); 512*09573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM,Mat); 513*09573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM,Mat); 514*09573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM,Mat); 515*09573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM,Mat); 516*09573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM,Mat); 517*09573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM,Mat); 518*09573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM,Mat); 519*09573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPISBAIJ(DM,Mat); 52047c6ae99SBarry Smith 52147c6ae99SBarry Smith #undef __FUNCT__ 52247c6ae99SBarry Smith #define __FUNCT__ "MatSetDA" 52347c6ae99SBarry Smith /*@ 524aa219208SBarry Smith MatSetDA - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix 52547c6ae99SBarry Smith 52647c6ae99SBarry Smith Logically Collective on Mat 52747c6ae99SBarry Smith 52847c6ae99SBarry Smith Input Parameters: 52947c6ae99SBarry Smith + mat - the matrix 53047c6ae99SBarry Smith - da - the da 53147c6ae99SBarry Smith 53247c6ae99SBarry Smith Level: intermediate 53347c6ae99SBarry Smith 53447c6ae99SBarry Smith @*/ 5359a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT MatSetDA(Mat mat,DM da) 53647c6ae99SBarry Smith { 53747c6ae99SBarry Smith PetscErrorCode ierr; 53847c6ae99SBarry Smith 53947c6ae99SBarry Smith PetscFunctionBegin; 54047c6ae99SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 54147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 5429a42bb27SBarry Smith ierr = PetscTryMethod(mat,"MatSetDA_C",(Mat,DM),(mat,da));CHKERRQ(ierr); 54347c6ae99SBarry Smith PetscFunctionReturn(0); 54447c6ae99SBarry Smith } 54547c6ae99SBarry Smith 54647c6ae99SBarry Smith EXTERN_C_BEGIN 54747c6ae99SBarry Smith #undef __FUNCT__ 54847c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA" 54947c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT MatView_MPI_DA(Mat A,PetscViewer viewer) 55047c6ae99SBarry Smith { 5519a42bb27SBarry Smith DM da; 55247c6ae99SBarry Smith PetscErrorCode ierr; 55347c6ae99SBarry Smith const char *prefix; 55447c6ae99SBarry Smith Mat Anatural; 55547c6ae99SBarry Smith AO ao; 55647c6ae99SBarry Smith PetscInt rstart,rend,*petsc,i; 55747c6ae99SBarry Smith IS is; 55847c6ae99SBarry Smith MPI_Comm comm; 55947c6ae99SBarry Smith 56047c6ae99SBarry Smith PetscFunctionBegin; 56147c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 562aa219208SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"DMDA",(PetscObject*)&da);CHKERRQ(ierr); 563aa219208SBarry Smith if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA"); 56447c6ae99SBarry Smith 565aa219208SBarry Smith ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr); 56647c6ae99SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 56747c6ae99SBarry Smith ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr); 56847c6ae99SBarry Smith for (i=rstart; i<rend; i++) petsc[i-rstart] = i; 56947c6ae99SBarry Smith ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr); 57047c6ae99SBarry Smith ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 57147c6ae99SBarry Smith 57247c6ae99SBarry Smith /* call viewer on natural ordering */ 57347c6ae99SBarry Smith ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr); 57447c6ae99SBarry Smith ierr = ISDestroy(is);CHKERRQ(ierr); 57547c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr); 57647c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr); 57747c6ae99SBarry Smith ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr); 57847c6ae99SBarry Smith ierr = MatView(Anatural,viewer);CHKERRQ(ierr); 57947c6ae99SBarry Smith ierr = MatDestroy(Anatural);CHKERRQ(ierr); 58047c6ae99SBarry Smith PetscFunctionReturn(0); 58147c6ae99SBarry Smith } 58247c6ae99SBarry Smith EXTERN_C_END 58347c6ae99SBarry Smith 58447c6ae99SBarry Smith EXTERN_C_BEGIN 58547c6ae99SBarry Smith #undef __FUNCT__ 58647c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA" 58747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT MatLoad_MPI_DA(Mat A,PetscViewer viewer) 58847c6ae99SBarry Smith { 5899a42bb27SBarry Smith DM da; 59047c6ae99SBarry Smith PetscErrorCode ierr; 59147c6ae99SBarry Smith Mat Anatural,Aapp; 59247c6ae99SBarry Smith AO ao; 59347c6ae99SBarry Smith PetscInt rstart,rend,*app,i; 59447c6ae99SBarry Smith IS is; 59547c6ae99SBarry Smith MPI_Comm comm; 59647c6ae99SBarry Smith 59747c6ae99SBarry Smith PetscFunctionBegin; 59847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 599aa219208SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"DMDA",(PetscObject*)&da);CHKERRQ(ierr); 600aa219208SBarry Smith if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA"); 60147c6ae99SBarry Smith 60247c6ae99SBarry Smith /* Load the matrix in natural ordering */ 60347c6ae99SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr); 60447c6ae99SBarry Smith ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr); 60547c6ae99SBarry Smith ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 60647c6ae99SBarry Smith ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr); 60747c6ae99SBarry Smith 60847c6ae99SBarry Smith /* Map natural ordering to application ordering and create IS */ 609aa219208SBarry Smith ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr); 61047c6ae99SBarry Smith ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr); 61147c6ae99SBarry Smith ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr); 61247c6ae99SBarry Smith for (i=rstart; i<rend; i++) app[i-rstart] = i; 61347c6ae99SBarry Smith ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr); 61447c6ae99SBarry Smith ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 61547c6ae99SBarry Smith 61647c6ae99SBarry Smith /* Do permutation and replace header */ 61747c6ae99SBarry Smith ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr); 61847c6ae99SBarry Smith ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr); 61947c6ae99SBarry Smith ierr = ISDestroy(is);CHKERRQ(ierr); 62047c6ae99SBarry Smith ierr = MatDestroy(Anatural);CHKERRQ(ierr); 62147c6ae99SBarry Smith PetscFunctionReturn(0); 62247c6ae99SBarry Smith } 62347c6ae99SBarry Smith EXTERN_C_END 62447c6ae99SBarry Smith 62547c6ae99SBarry Smith #undef __FUNCT__ 62694013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA" 62794013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix_DA(DM da, const MatType mtype,Mat *J) 62847c6ae99SBarry Smith { 62947c6ae99SBarry Smith PetscErrorCode ierr; 63047c6ae99SBarry Smith PetscInt dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P; 63147c6ae99SBarry Smith Mat A; 63247c6ae99SBarry Smith MPI_Comm comm; 63347c6ae99SBarry Smith const MatType Atype; 63447c6ae99SBarry Smith void (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL; 63547c6ae99SBarry Smith MatType ttype[256]; 63647c6ae99SBarry Smith PetscBool flg; 63747c6ae99SBarry Smith PetscMPIInt size; 63847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 63947c6ae99SBarry Smith 64047c6ae99SBarry Smith PetscFunctionBegin; 64147c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES 64247c6ae99SBarry Smith ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 64347c6ae99SBarry Smith #endif 6445da5aae0SJed Brown if (!mtype) mtype = MATAIJ; 64547c6ae99SBarry Smith ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr); 646aa219208SBarry Smith ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA options","Mat");CHKERRQ(ierr); 64747c6ae99SBarry Smith ierr = PetscOptionsList("-da_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr); 64847c6ae99SBarry Smith ierr = PetscOptionsEnd(); 64947c6ae99SBarry Smith 65047c6ae99SBarry Smith /* 65147c6ae99SBarry Smith m 65247c6ae99SBarry Smith ------------------------------------------------------ 65347c6ae99SBarry Smith | | 65447c6ae99SBarry Smith | | 65547c6ae99SBarry Smith | ---------------------- | 65647c6ae99SBarry Smith | | | | 65747c6ae99SBarry Smith n | ny | | | 65847c6ae99SBarry Smith | | | | 65947c6ae99SBarry Smith | .--------------------- | 66047c6ae99SBarry Smith | (xs,ys) nx | 66147c6ae99SBarry Smith | . | 66247c6ae99SBarry Smith | (gxs,gys) | 66347c6ae99SBarry Smith | | 66447c6ae99SBarry Smith ----------------------------------------------------- 66547c6ae99SBarry Smith */ 66647c6ae99SBarry Smith 66747c6ae99SBarry Smith /* 66847c6ae99SBarry Smith nc - number of components per grid point 66947c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 67047c6ae99SBarry Smith 67147c6ae99SBarry Smith */ 672aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0);CHKERRQ(ierr); 673aa219208SBarry Smith ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr); 67447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 67547c6ae99SBarry Smith ierr = MatCreate(comm,&A);CHKERRQ(ierr); 67647c6ae99SBarry Smith ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr); 67747c6ae99SBarry Smith ierr = MatSetType(A,(const MatType)ttype);CHKERRQ(ierr); 67847c6ae99SBarry Smith ierr = MatSetDA(A,da);CHKERRQ(ierr); 67947c6ae99SBarry Smith ierr = MatSetFromOptions(A);CHKERRQ(ierr); 68047c6ae99SBarry Smith ierr = MatGetType(A,&Atype);CHKERRQ(ierr); 68147c6ae99SBarry Smith /* 682aa219208SBarry Smith We do not provide a getmatrix function in the DMDA operations because 683aa219208SBarry Smith the basic DMDA does not know about matrices. We think of DMDA as being more 68447c6ae99SBarry Smith more low-level than matrices. This is kind of cheating but, cause sometimes 685aa219208SBarry Smith we think of DMDA has higher level than matrices. 68647c6ae99SBarry Smith 68747c6ae99SBarry Smith We could switch based on Atype (or mtype), but we do not since the 68847c6ae99SBarry Smith specialized setting routines depend only the particular preallocation 68947c6ae99SBarry Smith details of the matrix, not the type itself. 69047c6ae99SBarry Smith */ 69147c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 69247c6ae99SBarry Smith if (!aij) { 69347c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 69447c6ae99SBarry Smith } 69547c6ae99SBarry Smith if (!aij) { 69647c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr); 69747c6ae99SBarry Smith if (!baij) { 69847c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr); 69947c6ae99SBarry Smith } 70047c6ae99SBarry Smith if (!baij){ 70147c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr); 70247c6ae99SBarry Smith if (!sbaij) { 70347c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr); 70447c6ae99SBarry Smith } 70547c6ae99SBarry Smith } 70647c6ae99SBarry Smith } 70747c6ae99SBarry Smith if (aij) { 70847c6ae99SBarry Smith if (dim == 1) { 70994013140SBarry Smith ierr = DMGetMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr); 71047c6ae99SBarry Smith } else if (dim == 2) { 71147c6ae99SBarry Smith if (dd->ofill) { 71294013140SBarry Smith ierr = DMGetMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr); 71347c6ae99SBarry Smith } else { 71494013140SBarry Smith ierr = DMGetMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr); 71547c6ae99SBarry Smith } 71647c6ae99SBarry Smith } else if (dim == 3) { 71747c6ae99SBarry Smith if (dd->ofill) { 71894013140SBarry Smith ierr = DMGetMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr); 71947c6ae99SBarry Smith } else { 72094013140SBarry Smith ierr = DMGetMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr); 72147c6ae99SBarry Smith } 72247c6ae99SBarry Smith } 72347c6ae99SBarry Smith } else if (baij) { 72447c6ae99SBarry Smith if (dim == 2) { 72594013140SBarry Smith ierr = DMGetMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr); 72647c6ae99SBarry Smith } else if (dim == 3) { 72794013140SBarry Smith ierr = DMGetMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr); 72847c6ae99SBarry Smith } else { 729b17742caSSean Farley SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \ 730b17742caSSean Farley "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim); 73147c6ae99SBarry Smith } 73247c6ae99SBarry Smith } else if (sbaij) { 73347c6ae99SBarry Smith if (dim == 2) { 73494013140SBarry Smith ierr = DMGetMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr); 73547c6ae99SBarry Smith } else if (dim == 3) { 73694013140SBarry Smith ierr = DMGetMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr); 73747c6ae99SBarry Smith } else { 738b17742caSSean Farley SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \ 739b17742caSSean Farley "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim); 74047c6ae99SBarry Smith } 74147c6ae99SBarry Smith } 742aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr); 74347c6ae99SBarry Smith ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr); 744aa219208SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"DMDA",(PetscObject)da);CHKERRQ(ierr); 74547c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 74647c6ae99SBarry Smith if (size > 1) { 74747c6ae99SBarry Smith /* change viewer to display matrix in natural ordering */ 74847c6ae99SBarry Smith ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr); 74947c6ae99SBarry Smith ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr); 75047c6ae99SBarry Smith } 75147c6ae99SBarry Smith *J = A; 75247c6ae99SBarry Smith PetscFunctionReturn(0); 75347c6ae99SBarry Smith } 75447c6ae99SBarry Smith 75547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 75647c6ae99SBarry Smith #undef __FUNCT__ 75794013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIAIJ" 75894013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM da,Mat J) 75947c6ae99SBarry Smith { 76047c6ae99SBarry Smith PetscErrorCode ierr; 76147c6ae99SBarry 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; 76247c6ae99SBarry Smith PetscInt lstart,lend,pstart,pend,*dnz,*onz; 76347c6ae99SBarry Smith MPI_Comm comm; 76447c6ae99SBarry Smith PetscScalar *values; 765aa219208SBarry Smith DMDAPeriodicType wrap; 76647c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 767aa219208SBarry Smith DMDAStencilType st; 76847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 76947c6ae99SBarry Smith 77047c6ae99SBarry Smith PetscFunctionBegin; 77147c6ae99SBarry Smith /* 77247c6ae99SBarry Smith nc - number of components per grid point 77347c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 77447c6ae99SBarry Smith 77547c6ae99SBarry Smith */ 776aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr); 77747c6ae99SBarry Smith col = 2*s + 1; 778aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 779aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 78047c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 78147c6ae99SBarry Smith 78247c6ae99SBarry Smith ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr); 7831411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 7841411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 78547c6ae99SBarry Smith 78647c6ae99SBarry Smith /* determine the matrix preallocation information */ 78747c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr); 78847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 78947c6ae99SBarry Smith 790aa219208SBarry Smith pstart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 791aa219208SBarry Smith pend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 79247c6ae99SBarry Smith 79347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 79447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 79547c6ae99SBarry Smith 796aa219208SBarry Smith lstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 797aa219208SBarry Smith lend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 79847c6ae99SBarry Smith 79947c6ae99SBarry Smith cnt = 0; 80047c6ae99SBarry Smith for (k=0; k<nc; k++) { 80147c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 80247c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 803aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */ 80447c6ae99SBarry Smith cols[cnt++] = k + nc*(slot + gnx*l + p); 80547c6ae99SBarry Smith } 80647c6ae99SBarry Smith } 80747c6ae99SBarry Smith } 80847c6ae99SBarry Smith rows[k] = k + nc*(slot); 80947c6ae99SBarry Smith } 810784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 81147c6ae99SBarry Smith } 81247c6ae99SBarry Smith } 81347c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 81447c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 81547c6ae99SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 81647c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 81747c6ae99SBarry Smith 818784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 819784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 82047c6ae99SBarry Smith 82147c6ae99SBarry Smith /* 82247c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 82347c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 82447c6ae99SBarry Smith PETSc ordering. 82547c6ae99SBarry Smith */ 82647c6ae99SBarry Smith if (!dd->prealloc_only) { 82747c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 82847c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 82947c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 83047c6ae99SBarry Smith 831aa219208SBarry Smith pstart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 832aa219208SBarry Smith pend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 83347c6ae99SBarry Smith 83447c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 83547c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 83647c6ae99SBarry Smith 837aa219208SBarry Smith lstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 838aa219208SBarry Smith lend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 83947c6ae99SBarry Smith 84047c6ae99SBarry Smith cnt = 0; 84147c6ae99SBarry Smith for (k=0; k<nc; k++) { 84247c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 84347c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 844aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */ 84547c6ae99SBarry Smith cols[cnt++] = k + nc*(slot + gnx*l + p); 84647c6ae99SBarry Smith } 84747c6ae99SBarry Smith } 84847c6ae99SBarry Smith } 84947c6ae99SBarry Smith rows[k] = k + nc*(slot); 85047c6ae99SBarry Smith } 85147c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 85247c6ae99SBarry Smith } 85347c6ae99SBarry Smith } 85447c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 85547c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 85647c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 85747c6ae99SBarry Smith } 85847c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 85947c6ae99SBarry Smith PetscFunctionReturn(0); 86047c6ae99SBarry Smith } 86147c6ae99SBarry Smith 86247c6ae99SBarry Smith #undef __FUNCT__ 86394013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIAIJ_Fill" 86494013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J) 86547c6ae99SBarry Smith { 86647c6ae99SBarry Smith PetscErrorCode ierr; 86747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 86847c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,row,col,cnt,l,p; 86947c6ae99SBarry Smith PetscInt lstart,lend,pstart,pend,*dnz,*onz; 87047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 87147c6ae99SBarry Smith PetscInt ifill_col,*ofill = dd->ofill, *dfill = dd->dfill; 87247c6ae99SBarry Smith MPI_Comm comm; 87347c6ae99SBarry Smith PetscScalar *values; 874aa219208SBarry Smith DMDAPeriodicType wrap; 87547c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 876aa219208SBarry Smith DMDAStencilType st; 87747c6ae99SBarry Smith 87847c6ae99SBarry Smith PetscFunctionBegin; 87947c6ae99SBarry Smith /* 88047c6ae99SBarry Smith nc - number of components per grid point 88147c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 88247c6ae99SBarry Smith 88347c6ae99SBarry Smith */ 884aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr); 88547c6ae99SBarry Smith col = 2*s + 1; 886aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 887aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 88847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 88947c6ae99SBarry Smith 89047c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 8911411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 8921411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 89347c6ae99SBarry Smith 89447c6ae99SBarry Smith /* determine the matrix preallocation information */ 89547c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr); 89647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 89747c6ae99SBarry Smith 898aa219208SBarry Smith pstart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 899aa219208SBarry Smith pend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 90047c6ae99SBarry Smith 90147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 90247c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 90347c6ae99SBarry Smith 904aa219208SBarry Smith lstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 905aa219208SBarry Smith lend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 90647c6ae99SBarry Smith 90747c6ae99SBarry Smith for (k=0; k<nc; k++) { 90847c6ae99SBarry Smith cnt = 0; 90947c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 91047c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 91147c6ae99SBarry Smith if (l || p) { 912aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star */ 91347c6ae99SBarry Smith for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) 91447c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p); 91547c6ae99SBarry Smith } 91647c6ae99SBarry Smith } else { 91747c6ae99SBarry Smith if (dfill) { 91847c6ae99SBarry Smith for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) 91947c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p); 92047c6ae99SBarry Smith } else { 92147c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 92247c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + gnx*l + p); 92347c6ae99SBarry Smith } 92447c6ae99SBarry Smith } 92547c6ae99SBarry Smith } 92647c6ae99SBarry Smith } 92747c6ae99SBarry Smith row = k + nc*(slot); 928784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 92947c6ae99SBarry Smith } 93047c6ae99SBarry Smith } 93147c6ae99SBarry Smith } 93247c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 93347c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 93447c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 935784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 936784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 93747c6ae99SBarry Smith 93847c6ae99SBarry Smith /* 93947c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 94047c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 94147c6ae99SBarry Smith PETSc ordering. 94247c6ae99SBarry Smith */ 94347c6ae99SBarry Smith if (!dd->prealloc_only) { 94447c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 94547c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 94647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 94747c6ae99SBarry Smith 948aa219208SBarry Smith pstart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 949aa219208SBarry Smith pend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 95047c6ae99SBarry Smith 95147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 95247c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 95347c6ae99SBarry Smith 954aa219208SBarry Smith lstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 955aa219208SBarry Smith lend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 95647c6ae99SBarry Smith 95747c6ae99SBarry Smith for (k=0; k<nc; k++) { 95847c6ae99SBarry Smith cnt = 0; 95947c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 96047c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 96147c6ae99SBarry Smith if (l || p) { 962aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star */ 96347c6ae99SBarry Smith for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) 96447c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p); 96547c6ae99SBarry Smith } 96647c6ae99SBarry Smith } else { 96747c6ae99SBarry Smith if (dfill) { 96847c6ae99SBarry Smith for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) 96947c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p); 97047c6ae99SBarry Smith } else { 97147c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 97247c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + gnx*l + p); 97347c6ae99SBarry Smith } 97447c6ae99SBarry Smith } 97547c6ae99SBarry Smith } 97647c6ae99SBarry Smith } 97747c6ae99SBarry Smith row = k + nc*(slot); 97847c6ae99SBarry Smith ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 97947c6ae99SBarry Smith } 98047c6ae99SBarry Smith } 98147c6ae99SBarry Smith } 98247c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 98347c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 98447c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 98547c6ae99SBarry Smith } 98647c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 98747c6ae99SBarry Smith PetscFunctionReturn(0); 98847c6ae99SBarry Smith } 98947c6ae99SBarry Smith 99047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 99147c6ae99SBarry Smith 99247c6ae99SBarry Smith #undef __FUNCT__ 99394013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIAIJ" 99494013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM da,Mat J) 99547c6ae99SBarry Smith { 99647c6ae99SBarry Smith PetscErrorCode ierr; 99747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 99847c6ae99SBarry Smith PetscInt m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL; 99947c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 100047c6ae99SBarry Smith MPI_Comm comm; 100147c6ae99SBarry Smith PetscScalar *values; 1002aa219208SBarry Smith DMDAPeriodicType wrap; 100347c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 1004aa219208SBarry Smith DMDAStencilType st; 100547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 100647c6ae99SBarry Smith 100747c6ae99SBarry Smith PetscFunctionBegin; 100847c6ae99SBarry Smith /* 100947c6ae99SBarry Smith nc - number of components per grid point 101047c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 101147c6ae99SBarry Smith 101247c6ae99SBarry Smith */ 1013aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr); 101447c6ae99SBarry Smith col = 2*s + 1; 101547c6ae99SBarry Smith 1016aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1017aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 101847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 101947c6ae99SBarry Smith 102047c6ae99SBarry Smith ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr); 10211411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 10221411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 102347c6ae99SBarry Smith 102447c6ae99SBarry Smith /* determine the matrix preallocation information */ 102547c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr); 102647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1027aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1028aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 102947c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1030aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1031aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 103247c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1033aa219208SBarry Smith kstart = DMDAZPeriodic(wrap) ? -s : (PetscMax(-s,-k)); 1034aa219208SBarry Smith kend = DMDAZPeriodic(wrap) ? s : (PetscMin(s,p-k-1)); 103547c6ae99SBarry Smith 103647c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 103747c6ae99SBarry Smith 103847c6ae99SBarry Smith cnt = 0; 103947c6ae99SBarry Smith for (l=0; l<nc; l++) { 104047c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 104147c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 104247c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1043aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 104447c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk); 104547c6ae99SBarry Smith } 104647c6ae99SBarry Smith } 104747c6ae99SBarry Smith } 104847c6ae99SBarry Smith } 104947c6ae99SBarry Smith rows[l] = l + nc*(slot); 105047c6ae99SBarry Smith } 1051784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 105247c6ae99SBarry Smith } 105347c6ae99SBarry Smith } 105447c6ae99SBarry Smith } 105547c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 105647c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 105747c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 105847c6ae99SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 1059784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1060784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 106147c6ae99SBarry Smith 106247c6ae99SBarry Smith /* 106347c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 106447c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 106547c6ae99SBarry Smith PETSc ordering. 106647c6ae99SBarry Smith */ 106747c6ae99SBarry Smith if (!dd->prealloc_only) { 106847c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 106947c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 107047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1071aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1072aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 107347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1074aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1075aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 107647c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1077aa219208SBarry Smith kstart = DMDAZPeriodic(wrap) ? -s : (PetscMax(-s,-k)); 1078aa219208SBarry Smith kend = DMDAZPeriodic(wrap) ? s : (PetscMin(s,p-k-1)); 107947c6ae99SBarry Smith 108047c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 108147c6ae99SBarry Smith 108247c6ae99SBarry Smith cnt = 0; 108347c6ae99SBarry Smith for (l=0; l<nc; l++) { 108447c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 108547c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 108647c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1087aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 108847c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk); 108947c6ae99SBarry Smith } 109047c6ae99SBarry Smith } 109147c6ae99SBarry Smith } 109247c6ae99SBarry Smith } 109347c6ae99SBarry Smith rows[l] = l + nc*(slot); 109447c6ae99SBarry Smith } 109547c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 109647c6ae99SBarry Smith } 109747c6ae99SBarry Smith } 109847c6ae99SBarry Smith } 109947c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 110047c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 110147c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 110247c6ae99SBarry Smith } 110347c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 110447c6ae99SBarry Smith PetscFunctionReturn(0); 110547c6ae99SBarry Smith } 110647c6ae99SBarry Smith 110747c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 110847c6ae99SBarry Smith 110947c6ae99SBarry Smith #undef __FUNCT__ 111094013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_1d_MPIAIJ" 111194013140SBarry Smith PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM da,Mat J) 111247c6ae99SBarry Smith { 111347c6ae99SBarry Smith PetscErrorCode ierr; 111447c6ae99SBarry Smith PetscInt xs,nx,i,i1,slot,gxs,gnx; 111547c6ae99SBarry Smith PetscInt m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l; 111647c6ae99SBarry Smith PetscInt istart,iend; 111747c6ae99SBarry Smith PetscScalar *values; 1118aa219208SBarry Smith DMDAPeriodicType wrap; 111947c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 112047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 112147c6ae99SBarry Smith 112247c6ae99SBarry Smith PetscFunctionBegin; 112347c6ae99SBarry Smith /* 112447c6ae99SBarry Smith nc - number of components per grid point 112547c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 112647c6ae99SBarry Smith 112747c6ae99SBarry Smith */ 1128aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&wrap,0);CHKERRQ(ierr); 112947c6ae99SBarry Smith col = 2*s + 1; 113047c6ae99SBarry Smith 1131aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr); 1132aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr); 113347c6ae99SBarry Smith 113447c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr); 113547c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr); 113647c6ae99SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 113747c6ae99SBarry Smith ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr); 113847c6ae99SBarry Smith 11391411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 11401411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 1141784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1142784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 114347c6ae99SBarry Smith 114447c6ae99SBarry Smith /* 114547c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 114647c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 114747c6ae99SBarry Smith PETSc ordering. 114847c6ae99SBarry Smith */ 114947c6ae99SBarry Smith if (!dd->prealloc_only) { 115047c6ae99SBarry Smith ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 115147c6ae99SBarry Smith ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 115247c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 115347c6ae99SBarry Smith istart = PetscMax(-s,gxs - i); 115447c6ae99SBarry Smith iend = PetscMin(s,gxs + gnx - i - 1); 115547c6ae99SBarry Smith slot = i - gxs; 115647c6ae99SBarry Smith 115747c6ae99SBarry Smith cnt = 0; 115847c6ae99SBarry Smith for (l=0; l<nc; l++) { 115947c6ae99SBarry Smith for (i1=istart; i1<iend+1; i1++) { 116047c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + i1); 116147c6ae99SBarry Smith } 116247c6ae99SBarry Smith rows[l] = l + nc*(slot); 116347c6ae99SBarry Smith } 116447c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 116547c6ae99SBarry Smith } 116647c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 116747c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 116847c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 116947c6ae99SBarry Smith } 117047c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 117147c6ae99SBarry Smith PetscFunctionReturn(0); 117247c6ae99SBarry Smith } 117347c6ae99SBarry Smith 117447c6ae99SBarry Smith #undef __FUNCT__ 117594013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIBAIJ" 117694013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM da,Mat J) 117747c6ae99SBarry Smith { 117847c6ae99SBarry Smith PetscErrorCode ierr; 117947c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 118047c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz; 118147c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,ii,jj; 118247c6ae99SBarry Smith MPI_Comm comm; 118347c6ae99SBarry Smith PetscScalar *values; 1184aa219208SBarry Smith DMDAPeriodicType wrap; 1185aa219208SBarry Smith DMDAStencilType st; 118647c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 118747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 118847c6ae99SBarry Smith 118947c6ae99SBarry Smith PetscFunctionBegin; 119047c6ae99SBarry Smith /* 119147c6ae99SBarry Smith nc - number of components per grid point 119247c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 119347c6ae99SBarry Smith */ 1194aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr); 119547c6ae99SBarry Smith col = 2*s + 1; 119647c6ae99SBarry Smith 1197aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 1198aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 119947c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 120047c6ae99SBarry Smith 120147c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 120247c6ae99SBarry Smith 12031411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 12041411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 120547c6ae99SBarry Smith 120647c6ae99SBarry Smith /* determine the matrix preallocation information */ 120747c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr); 120847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1209aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1210aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 121147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1212aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1213aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 121447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 121547c6ae99SBarry Smith 121647c6ae99SBarry Smith /* Find block columns in block row */ 121747c6ae99SBarry Smith cnt = 0; 121847c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 121947c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1220aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */ 122147c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 122247c6ae99SBarry Smith } 122347c6ae99SBarry Smith } 122447c6ae99SBarry Smith } 1225784ac674SJed Brown ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr); 122647c6ae99SBarry Smith } 122747c6ae99SBarry Smith } 122847c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 122947c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 123047c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 123147c6ae99SBarry Smith 1232784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1233784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 123447c6ae99SBarry Smith 123547c6ae99SBarry Smith /* 123647c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 123747c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 123847c6ae99SBarry Smith PETSc ordering. 123947c6ae99SBarry Smith */ 124047c6ae99SBarry Smith if (!dd->prealloc_only) { 124147c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 124247c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 124347c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1244aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1245aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 124647c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1247aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1248aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 124947c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 125047c6ae99SBarry Smith cnt = 0; 125147c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 125247c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1253aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */ 125447c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 125547c6ae99SBarry Smith } 125647c6ae99SBarry Smith } 125747c6ae99SBarry Smith } 125847c6ae99SBarry Smith ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 125947c6ae99SBarry Smith } 126047c6ae99SBarry Smith } 126147c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 126247c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126347c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126447c6ae99SBarry Smith } 126547c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 126647c6ae99SBarry Smith PetscFunctionReturn(0); 126747c6ae99SBarry Smith } 126847c6ae99SBarry Smith 126947c6ae99SBarry Smith #undef __FUNCT__ 127094013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIBAIJ" 127194013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM da,Mat J) 127247c6ae99SBarry Smith { 127347c6ae99SBarry Smith PetscErrorCode ierr; 127447c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 127547c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz; 127647c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 127747c6ae99SBarry Smith MPI_Comm comm; 127847c6ae99SBarry Smith PetscScalar *values; 1279aa219208SBarry Smith DMDAPeriodicType wrap; 1280aa219208SBarry Smith DMDAStencilType st; 128147c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 128247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 128347c6ae99SBarry Smith 128447c6ae99SBarry Smith PetscFunctionBegin; 128547c6ae99SBarry Smith /* 128647c6ae99SBarry Smith nc - number of components per grid point 128747c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 128847c6ae99SBarry Smith 128947c6ae99SBarry Smith */ 1290aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr); 129147c6ae99SBarry Smith col = 2*s + 1; 129247c6ae99SBarry Smith 1293aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1294aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 129547c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 129647c6ae99SBarry Smith 129747c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr); 129847c6ae99SBarry Smith 12991411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 13001411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 130147c6ae99SBarry Smith 130247c6ae99SBarry Smith /* determine the matrix preallocation information */ 130347c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr); 130447c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1305aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1306aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 130747c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1308aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1309aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 131047c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1311aa219208SBarry Smith kstart = DMDAZPeriodic(wrap) ? -s : (PetscMax(-s,-k)); 1312aa219208SBarry Smith kend = DMDAZPeriodic(wrap) ? s : (PetscMin(s,p-k-1)); 131347c6ae99SBarry Smith 131447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 131547c6ae99SBarry Smith 131647c6ae99SBarry Smith /* Find block columns in block row */ 131747c6ae99SBarry Smith cnt = 0; 131847c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 131947c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 132047c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1321aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 132247c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 132347c6ae99SBarry Smith } 132447c6ae99SBarry Smith } 132547c6ae99SBarry Smith } 132647c6ae99SBarry Smith } 1327784ac674SJed Brown ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr); 132847c6ae99SBarry Smith } 132947c6ae99SBarry Smith } 133047c6ae99SBarry Smith } 133147c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 133247c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 133347c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 133447c6ae99SBarry Smith 1335784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1336784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 133747c6ae99SBarry Smith 133847c6ae99SBarry Smith /* 133947c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 134047c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 134147c6ae99SBarry Smith PETSc ordering. 134247c6ae99SBarry Smith */ 134347c6ae99SBarry Smith if (!dd->prealloc_only) { 134447c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 134547c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 134647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1347aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1348aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 134947c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1350aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1351aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 135247c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1353aa219208SBarry Smith kstart = DMDAZPeriodic(wrap) ? -s : (PetscMax(-s,-k)); 1354aa219208SBarry Smith kend = DMDAZPeriodic(wrap) ? s : (PetscMin(s,p-k-1)); 135547c6ae99SBarry Smith 135647c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 135747c6ae99SBarry Smith 135847c6ae99SBarry Smith cnt = 0; 135947c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 136047c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 136147c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1362aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 136347c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 136447c6ae99SBarry Smith } 136547c6ae99SBarry Smith } 136647c6ae99SBarry Smith } 136747c6ae99SBarry Smith } 136847c6ae99SBarry Smith ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 136947c6ae99SBarry Smith } 137047c6ae99SBarry Smith } 137147c6ae99SBarry Smith } 137247c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 137347c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 137447c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 137547c6ae99SBarry Smith } 137647c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 137747c6ae99SBarry Smith PetscFunctionReturn(0); 137847c6ae99SBarry Smith } 137947c6ae99SBarry Smith 138047c6ae99SBarry Smith #undef __FUNCT__ 138147c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular" 138247c6ae99SBarry Smith /* 138347c6ae99SBarry Smith This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to 138447c6ae99SBarry Smith identify in the local ordering with periodic domain. 138547c6ae99SBarry Smith */ 138647c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[]) 138747c6ae99SBarry Smith { 138847c6ae99SBarry Smith PetscErrorCode ierr; 138947c6ae99SBarry Smith PetscInt i,n; 139047c6ae99SBarry Smith 139147c6ae99SBarry Smith PetscFunctionBegin; 139247c6ae99SBarry Smith ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr); 139347c6ae99SBarry Smith ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr); 139447c6ae99SBarry Smith for (i=0,n=0; i<*cnt; i++) { 139547c6ae99SBarry Smith if (col[i] >= *row) col[n++] = col[i]; 139647c6ae99SBarry Smith } 139747c6ae99SBarry Smith *cnt = n; 139847c6ae99SBarry Smith PetscFunctionReturn(0); 139947c6ae99SBarry Smith } 140047c6ae99SBarry Smith 140147c6ae99SBarry Smith #undef __FUNCT__ 140294013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPISBAIJ" 140394013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM da,Mat J) 140447c6ae99SBarry Smith { 140547c6ae99SBarry Smith PetscErrorCode ierr; 140647c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 140747c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz; 140847c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,ii,jj; 140947c6ae99SBarry Smith MPI_Comm comm; 141047c6ae99SBarry Smith PetscScalar *values; 1411aa219208SBarry Smith DMDAPeriodicType wrap; 1412aa219208SBarry Smith DMDAStencilType st; 141347c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 141447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 141547c6ae99SBarry Smith 141647c6ae99SBarry Smith PetscFunctionBegin; 141747c6ae99SBarry Smith /* 141847c6ae99SBarry Smith nc - number of components per grid point 141947c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 142047c6ae99SBarry Smith */ 1421aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr); 142247c6ae99SBarry Smith col = 2*s + 1; 142347c6ae99SBarry Smith 1424aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 1425aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 142647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 142747c6ae99SBarry Smith 142847c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 142947c6ae99SBarry Smith 14301411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 14311411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 143247c6ae99SBarry Smith 143347c6ae99SBarry Smith /* determine the matrix preallocation information */ 143447c6ae99SBarry Smith ierr = MatPreallocateSymmetricInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr); 143547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1436aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1437aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 143847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1439aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1440aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 144147c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 144247c6ae99SBarry Smith 144347c6ae99SBarry Smith /* Find block columns in block row */ 144447c6ae99SBarry Smith cnt = 0; 144547c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 144647c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1447aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { 144847c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 144947c6ae99SBarry Smith } 145047c6ae99SBarry Smith } 145147c6ae99SBarry Smith } 145247c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 145347c6ae99SBarry Smith ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr); 145447c6ae99SBarry Smith } 145547c6ae99SBarry Smith } 145647c6ae99SBarry Smith ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 145747c6ae99SBarry Smith ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 145847c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 145947c6ae99SBarry Smith 1460784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1461784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 146247c6ae99SBarry Smith 146347c6ae99SBarry Smith /* 146447c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 146547c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 146647c6ae99SBarry Smith PETSc ordering. 146747c6ae99SBarry Smith */ 146847c6ae99SBarry Smith if (!dd->prealloc_only) { 146947c6ae99SBarry Smith ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 147047c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 147147c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1472aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1473aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 147447c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1475aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1476aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 147747c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 147847c6ae99SBarry Smith 147947c6ae99SBarry Smith /* Find block columns in block row */ 148047c6ae99SBarry Smith cnt = 0; 148147c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 148247c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1483aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { 148447c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 148547c6ae99SBarry Smith } 148647c6ae99SBarry Smith } 148747c6ae99SBarry Smith } 148847c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 148947c6ae99SBarry Smith ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 149047c6ae99SBarry Smith } 149147c6ae99SBarry Smith } 149247c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 149347c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 149447c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 149547c6ae99SBarry Smith } 149647c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 149747c6ae99SBarry Smith PetscFunctionReturn(0); 149847c6ae99SBarry Smith } 149947c6ae99SBarry Smith 150047c6ae99SBarry Smith #undef __FUNCT__ 150194013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPISBAIJ" 150294013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPISBAIJ(DM da,Mat J) 150347c6ae99SBarry Smith { 150447c6ae99SBarry Smith PetscErrorCode ierr; 150547c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 150647c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz; 150747c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 150847c6ae99SBarry Smith MPI_Comm comm; 150947c6ae99SBarry Smith PetscScalar *values; 1510aa219208SBarry Smith DMDAPeriodicType wrap; 1511aa219208SBarry Smith DMDAStencilType st; 151247c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 151347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 151447c6ae99SBarry Smith 151547c6ae99SBarry Smith PetscFunctionBegin; 151647c6ae99SBarry Smith /* 151747c6ae99SBarry Smith nc - number of components per grid point 151847c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 151947c6ae99SBarry Smith */ 1520aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr); 152147c6ae99SBarry Smith col = 2*s + 1; 152247c6ae99SBarry Smith 1523aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1524aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 152547c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 152647c6ae99SBarry Smith 152747c6ae99SBarry Smith /* create the matrix */ 152847c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr); 152947c6ae99SBarry Smith 15301411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 15311411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 153247c6ae99SBarry Smith 153347c6ae99SBarry Smith /* determine the matrix preallocation information */ 153447c6ae99SBarry Smith ierr = MatPreallocateSymmetricInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr); 153547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1536aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1537aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 153847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1539aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1540aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 154147c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1542aa219208SBarry Smith kstart = DMDAZPeriodic(wrap) ? -s : (PetscMax(-s,-k)); 1543aa219208SBarry Smith kend = DMDAZPeriodic(wrap) ? s : (PetscMin(s,p-k-1)); 154447c6ae99SBarry Smith 154547c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 154647c6ae99SBarry Smith 154747c6ae99SBarry Smith /* Find block columns in block row */ 154847c6ae99SBarry Smith cnt = 0; 154947c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 155047c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 155147c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1552aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) { 155347c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 155447c6ae99SBarry Smith } 155547c6ae99SBarry Smith } 155647c6ae99SBarry Smith } 155747c6ae99SBarry Smith } 155847c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 155947c6ae99SBarry Smith ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr); 156047c6ae99SBarry Smith } 156147c6ae99SBarry Smith } 156247c6ae99SBarry Smith } 156347c6ae99SBarry Smith ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 156447c6ae99SBarry Smith ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 156547c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 156647c6ae99SBarry Smith 1567784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1568784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 156947c6ae99SBarry Smith 157047c6ae99SBarry Smith /* 157147c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 157247c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 157347c6ae99SBarry Smith PETSc ordering. 157447c6ae99SBarry Smith */ 157547c6ae99SBarry Smith if (!dd->prealloc_only) { 157647c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 157747c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 157847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1579aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1580aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 158147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1582aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1583aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 158447c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1585aa219208SBarry Smith kstart = DMDAZPeriodic(wrap) ? -s : (PetscMax(-s,-k)); 1586aa219208SBarry Smith kend = DMDAZPeriodic(wrap) ? s : (PetscMin(s,p-k-1)); 158747c6ae99SBarry Smith 158847c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 158947c6ae99SBarry Smith 159047c6ae99SBarry Smith cnt = 0; 159147c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 159247c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 159347c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1594aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) { 159547c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 159647c6ae99SBarry Smith } 159747c6ae99SBarry Smith } 159847c6ae99SBarry Smith } 159947c6ae99SBarry Smith } 160047c6ae99SBarry Smith ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr); 160147c6ae99SBarry Smith ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 160247c6ae99SBarry Smith } 160347c6ae99SBarry Smith } 160447c6ae99SBarry Smith } 160547c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 160647c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 160747c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 160847c6ae99SBarry Smith } 160947c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 161047c6ae99SBarry Smith PetscFunctionReturn(0); 161147c6ae99SBarry Smith } 161247c6ae99SBarry Smith 161347c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 161447c6ae99SBarry Smith 161547c6ae99SBarry Smith #undef __FUNCT__ 161694013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIAIJ_Fill" 161794013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J) 161847c6ae99SBarry Smith { 161947c6ae99SBarry Smith PetscErrorCode ierr; 162047c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 162147c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz; 162247c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 162347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 162447c6ae99SBarry Smith PetscInt ifill_col,*dfill = dd->dfill,*ofill = dd->ofill; 162547c6ae99SBarry Smith MPI_Comm comm; 162647c6ae99SBarry Smith PetscScalar *values; 1627aa219208SBarry Smith DMDAPeriodicType wrap; 162847c6ae99SBarry Smith ISLocalToGlobalMapping ltog,ltogb; 1629aa219208SBarry Smith DMDAStencilType st; 163047c6ae99SBarry Smith 163147c6ae99SBarry Smith PetscFunctionBegin; 163247c6ae99SBarry Smith /* 163347c6ae99SBarry Smith nc - number of components per grid point 163447c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 163547c6ae99SBarry Smith 163647c6ae99SBarry Smith */ 1637aa219208SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr); 163847c6ae99SBarry Smith col = 2*s + 1; 1639aa219208SBarry Smith if (DMDAXPeriodic(wrap) && (m % col)){ 164047c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\ 164147c6ae99SBarry Smith by 2*stencil_width + 1\n"); 164247c6ae99SBarry Smith } 1643aa219208SBarry Smith if (DMDAYPeriodic(wrap) && (n % col)){ 164447c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\ 164547c6ae99SBarry Smith by 2*stencil_width + 1\n"); 164647c6ae99SBarry Smith } 1647aa219208SBarry Smith if (DMDAZPeriodic(wrap) && (p % col)){ 164847c6ae99SBarry Smith SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\ 164947c6ae99SBarry Smith by 2*stencil_width + 1\n"); 165047c6ae99SBarry Smith } 165147c6ae99SBarry Smith 1652aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1653aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 165447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 165547c6ae99SBarry Smith 165647c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr); 16571411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 16581411c6eeSJed Brown ierr = DMGetLocalToGlobalMappingBlock(da,<ogb);CHKERRQ(ierr); 165947c6ae99SBarry Smith 166047c6ae99SBarry Smith /* determine the matrix preallocation information */ 166147c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr); 166247c6ae99SBarry Smith 166347c6ae99SBarry Smith 166447c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1665aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1666aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 166747c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1668aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1669aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 167047c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1671aa219208SBarry Smith kstart = DMDAZPeriodic(wrap) ? -s : (PetscMax(-s,-k)); 1672aa219208SBarry Smith kend = DMDAZPeriodic(wrap) ? s : (PetscMin(s,p-k-1)); 167347c6ae99SBarry Smith 167447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 167547c6ae99SBarry Smith 167647c6ae99SBarry Smith for (l=0; l<nc; l++) { 167747c6ae99SBarry Smith cnt = 0; 167847c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 167947c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 168047c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 168147c6ae99SBarry Smith if (ii || jj || kk) { 1682aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 168347c6ae99SBarry Smith for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) 168447c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 168547c6ae99SBarry Smith } 168647c6ae99SBarry Smith } else { 168747c6ae99SBarry Smith if (dfill) { 168847c6ae99SBarry Smith for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) 168947c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 169047c6ae99SBarry Smith } else { 169147c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 169247c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk); 169347c6ae99SBarry Smith } 169447c6ae99SBarry Smith } 169547c6ae99SBarry Smith } 169647c6ae99SBarry Smith } 169747c6ae99SBarry Smith } 169847c6ae99SBarry Smith row = l + nc*(slot); 1699784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 170047c6ae99SBarry Smith } 170147c6ae99SBarry Smith } 170247c6ae99SBarry Smith } 170347c6ae99SBarry Smith } 170447c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 170547c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 170647c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1707784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1708784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr); 170947c6ae99SBarry Smith 171047c6ae99SBarry Smith /* 171147c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 171247c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 171347c6ae99SBarry Smith PETSc ordering. 171447c6ae99SBarry Smith */ 171547c6ae99SBarry Smith if (!dd->prealloc_only) { 171647c6ae99SBarry Smith ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr); 171747c6ae99SBarry Smith ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr); 171847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1719aa219208SBarry Smith istart = DMDAXPeriodic(wrap) ? -s : (PetscMax(-s,-i)); 1720aa219208SBarry Smith iend = DMDAXPeriodic(wrap) ? s : (PetscMin(s,m-i-1)); 172147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1722aa219208SBarry Smith jstart = DMDAYPeriodic(wrap) ? -s : (PetscMax(-s,-j)); 1723aa219208SBarry Smith jend = DMDAYPeriodic(wrap) ? s : (PetscMin(s,n-j-1)); 172447c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1725aa219208SBarry Smith kstart = DMDAZPeriodic(wrap) ? -s : (PetscMax(-s,-k)); 1726aa219208SBarry Smith kend = DMDAZPeriodic(wrap) ? s : (PetscMin(s,p-k-1)); 172747c6ae99SBarry Smith 172847c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 172947c6ae99SBarry Smith 173047c6ae99SBarry Smith for (l=0; l<nc; l++) { 173147c6ae99SBarry Smith cnt = 0; 173247c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 173347c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 173447c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 173547c6ae99SBarry Smith if (ii || jj || kk) { 1736aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 173747c6ae99SBarry Smith for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) 173847c6ae99SBarry Smith cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 173947c6ae99SBarry Smith } 174047c6ae99SBarry Smith } else { 174147c6ae99SBarry Smith if (dfill) { 174247c6ae99SBarry Smith for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) 174347c6ae99SBarry Smith cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 174447c6ae99SBarry Smith } else { 174547c6ae99SBarry Smith for (ifill_col=0; ifill_col<nc; ifill_col++) 174647c6ae99SBarry Smith cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk); 174747c6ae99SBarry Smith } 174847c6ae99SBarry Smith } 174947c6ae99SBarry Smith } 175047c6ae99SBarry Smith } 175147c6ae99SBarry Smith } 175247c6ae99SBarry Smith row = l + nc*(slot); 175347c6ae99SBarry Smith ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 175447c6ae99SBarry Smith } 175547c6ae99SBarry Smith } 175647c6ae99SBarry Smith } 175747c6ae99SBarry Smith } 175847c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 175947c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 176047c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 176147c6ae99SBarry Smith } 176247c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 176347c6ae99SBarry Smith PetscFunctionReturn(0); 176447c6ae99SBarry Smith } 1765