147c6ae99SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 307475bc1SBarry Smith #include <petscmat.h> 447c6ae99SBarry Smith 5e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring*); 6e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring*); 7e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring*); 8e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM,ISColoringType,ISColoring*); 947c6ae99SBarry Smith 1047c6ae99SBarry Smith /* 1147c6ae99SBarry Smith For ghost i that may be negative or greater than the upper bound this 1247c6ae99SBarry Smith maps it into the 0:m-1 range using periodicity 1347c6ae99SBarry Smith */ 1447c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i : ((i >= m) ? i-m : i)) 1547c6ae99SBarry Smith 16ce308e1dSBarry Smith static PetscErrorCode DMDASetBlockFills_Private(const PetscInt *dfill,PetscInt w,PetscInt **rfill) 1747c6ae99SBarry Smith { 1847c6ae99SBarry Smith PetscErrorCode ierr; 1947c6ae99SBarry Smith PetscInt i,j,nz,*fill; 2047c6ae99SBarry Smith 2147c6ae99SBarry Smith PetscFunctionBegin; 2247c6ae99SBarry Smith if (!dfill) PetscFunctionReturn(0); 2347c6ae99SBarry Smith 2447c6ae99SBarry Smith /* count number nonzeros */ 2547c6ae99SBarry Smith nz = 0; 2647c6ae99SBarry Smith for (i=0; i<w; i++) { 2747c6ae99SBarry Smith for (j=0; j<w; j++) { 2847c6ae99SBarry Smith if (dfill[w*i+j]) nz++; 2947c6ae99SBarry Smith } 3047c6ae99SBarry Smith } 31854ce69bSBarry Smith ierr = PetscMalloc1(nz + w + 1,&fill);CHKERRQ(ierr); 3247c6ae99SBarry Smith /* construct modified CSR storage of nonzero structure */ 33ce308e1dSBarry Smith /* fill[0 -- w] marks starts of each row of column indices (and end of last row) 34ce308e1dSBarry Smith so fill[1] - fill[0] gives number of nonzeros in first row etc */ 3547c6ae99SBarry Smith nz = w + 1; 3647c6ae99SBarry Smith for (i=0; i<w; i++) { 3747c6ae99SBarry Smith fill[i] = nz; 3847c6ae99SBarry Smith for (j=0; j<w; j++) { 3947c6ae99SBarry Smith if (dfill[w*i+j]) { 4047c6ae99SBarry Smith fill[nz] = j; 4147c6ae99SBarry Smith nz++; 4247c6ae99SBarry Smith } 4347c6ae99SBarry Smith } 4447c6ae99SBarry Smith } 4547c6ae99SBarry Smith fill[w] = nz; 4647c6ae99SBarry Smith 4747c6ae99SBarry Smith *rfill = fill; 4847c6ae99SBarry Smith PetscFunctionReturn(0); 4947c6ae99SBarry Smith } 5047c6ae99SBarry Smith 5147c6ae99SBarry Smith /*@ 52aa219208SBarry Smith DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem 53950540a4SJed Brown of the matrix returned by DMCreateMatrix(). 5447c6ae99SBarry Smith 55aa219208SBarry Smith Logically Collective on DMDA 5647c6ae99SBarry Smith 5747c6ae99SBarry Smith Input Parameter: 5847c6ae99SBarry Smith + da - the distributed array 590298fd71SBarry Smith . dfill - the fill pattern in the diagonal block (may be NULL, means use dense block) 6047c6ae99SBarry Smith - ofill - the fill pattern in the off-diagonal blocks 6147c6ae99SBarry Smith 6247c6ae99SBarry Smith 6347c6ae99SBarry Smith Level: developer 6447c6ae99SBarry Smith 6547c6ae99SBarry Smith Notes: This only makes sense when you are doing multicomponent problems but using the 6647c6ae99SBarry Smith MPIAIJ matrix format 6747c6ae99SBarry Smith 6847c6ae99SBarry Smith The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries 6947c6ae99SBarry Smith representing coupling and 0 entries for missing coupling. For example 7047c6ae99SBarry Smith $ dfill[9] = {1, 0, 0, 7147c6ae99SBarry Smith $ 1, 1, 0, 7247c6ae99SBarry Smith $ 0, 1, 1} 7347c6ae99SBarry Smith means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with 7447c6ae99SBarry Smith itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the 7547c6ae99SBarry Smith diagonal block). 7647c6ae99SBarry Smith 77aa219208SBarry Smith DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then 7847c6ae99SBarry Smith can be represented in the dfill, ofill format 7947c6ae99SBarry Smith 8047c6ae99SBarry Smith Contributed by Glenn Hammond 8147c6ae99SBarry Smith 828ddb5d8bSBarry Smith .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly() 8347c6ae99SBarry Smith 8447c6ae99SBarry Smith @*/ 85ce308e1dSBarry Smith PetscErrorCode DMDASetBlockFills(DM da,const PetscInt *dfill,const PetscInt *ofill) 8647c6ae99SBarry Smith { 8747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 8847c6ae99SBarry Smith PetscErrorCode ierr; 89ae4f298aSBarry Smith PetscInt i,k,cnt = 1; 9047c6ae99SBarry Smith 9147c6ae99SBarry Smith PetscFunctionBegin; 92aa219208SBarry Smith ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr); 93aa219208SBarry Smith ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr); 94ae4f298aSBarry Smith 95ae4f298aSBarry Smith /* ofillcount tracks the columns of ofill that have any nonzero in thems; the value in each location is the number of 96ae4f298aSBarry Smith columns to the left with any nonzeros in them plus 1 */ 971795a4d1SJed Brown ierr = PetscCalloc1(dd->w,&dd->ofillcols);CHKERRQ(ierr); 98ae4f298aSBarry Smith for (i=0; i<dd->w; i++) { 99ae4f298aSBarry Smith for (k=dd->ofill[i]; k<dd->ofill[i+1]; k++) dd->ofillcols[dd->ofill[k]] = 1; 100ae4f298aSBarry Smith } 101ae4f298aSBarry Smith for (i=0; i<dd->w; i++) { 102ae4f298aSBarry Smith if (dd->ofillcols[i]) { 103ae4f298aSBarry Smith dd->ofillcols[i] = cnt++; 104ae4f298aSBarry Smith } 105ae4f298aSBarry Smith } 10647c6ae99SBarry Smith PetscFunctionReturn(0); 10747c6ae99SBarry Smith } 10847c6ae99SBarry Smith 10947c6ae99SBarry Smith 110b412c318SBarry Smith PetscErrorCode DMCreateColoring_DA(DM da,ISColoringType ctype,ISColoring *coloring) 11147c6ae99SBarry Smith { 11247c6ae99SBarry Smith PetscErrorCode ierr; 11347c6ae99SBarry Smith PetscInt dim,m,n,p,nc; 114bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by,bz; 11547c6ae99SBarry Smith MPI_Comm comm; 11647c6ae99SBarry Smith PetscMPIInt size; 11747c6ae99SBarry Smith PetscBool isBAIJ; 11847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 11947c6ae99SBarry Smith 12047c6ae99SBarry Smith PetscFunctionBegin; 12147c6ae99SBarry Smith /* 12247c6ae99SBarry Smith m 12347c6ae99SBarry Smith ------------------------------------------------------ 12447c6ae99SBarry Smith | | 12547c6ae99SBarry Smith | | 12647c6ae99SBarry Smith | ---------------------- | 12747c6ae99SBarry Smith | | | | 12847c6ae99SBarry Smith n | yn | | | 12947c6ae99SBarry Smith | | | | 13047c6ae99SBarry Smith | .--------------------- | 13147c6ae99SBarry Smith | (xs,ys) xn | 13247c6ae99SBarry Smith | . | 13347c6ae99SBarry Smith | (gxs,gys) | 13447c6ae99SBarry Smith | | 13547c6ae99SBarry Smith ----------------------------------------------------- 13647c6ae99SBarry Smith */ 13747c6ae99SBarry Smith 13847c6ae99SBarry Smith /* 13947c6ae99SBarry Smith nc - number of components per grid point 14047c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 14147c6ae99SBarry Smith 14247c6ae99SBarry Smith */ 1431321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);CHKERRQ(ierr); 14447c6ae99SBarry Smith 14547c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 14647c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1475bdb020cSBarry Smith if (ctype == IS_COLORING_LOCAL) { 14847c6ae99SBarry Smith if (size == 1) { 14947c6ae99SBarry Smith ctype = IS_COLORING_GLOBAL; 15047c6ae99SBarry Smith } else if (dim > 1) { 151bff4a2f0SMatthew G. Knepley if ((m==1 && bx == DM_BOUNDARY_PERIODIC) || (n==1 && by == DM_BOUNDARY_PERIODIC) || (p==1 && bz == DM_BOUNDARY_PERIODIC)) { 1525bdb020cSBarry Smith SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"IS_COLORING_LOCAL cannot be used for periodic boundary condition having both ends of the domain on the same process"); 15347c6ae99SBarry Smith } 15447c6ae99SBarry Smith } 15547c6ae99SBarry Smith } 15647c6ae99SBarry Smith 157aa219208SBarry Smith /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ 15847c6ae99SBarry Smith matrices is for the blocks, not the individual matrix elements */ 159b412c318SBarry Smith ierr = PetscStrcmp(da->mattype,MATBAIJ,&isBAIJ);CHKERRQ(ierr); 160b412c318SBarry Smith if (!isBAIJ) {ierr = PetscStrcmp(da->mattype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);} 161b412c318SBarry Smith if (!isBAIJ) {ierr = PetscStrcmp(da->mattype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);} 16247c6ae99SBarry Smith if (isBAIJ) { 16347c6ae99SBarry Smith dd->w = 1; 16447c6ae99SBarry Smith dd->xs = dd->xs/nc; 16547c6ae99SBarry Smith dd->xe = dd->xe/nc; 16647c6ae99SBarry Smith dd->Xs = dd->Xs/nc; 16747c6ae99SBarry Smith dd->Xe = dd->Xe/nc; 16847c6ae99SBarry Smith } 16947c6ae99SBarry Smith 17047c6ae99SBarry Smith /* 171aa219208SBarry Smith We do not provide a getcoloring function in the DMDA operations because 172aa219208SBarry Smith the basic DMDA does not know about matrices. We think of DMDA as being more 17347c6ae99SBarry Smith more low-level then matrices. 17447c6ae99SBarry Smith */ 17547c6ae99SBarry Smith if (dim == 1) { 176e727c939SJed Brown ierr = DMCreateColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 17747c6ae99SBarry Smith } else if (dim == 2) { 178e727c939SJed Brown ierr = DMCreateColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 17947c6ae99SBarry Smith } else if (dim == 3) { 180e727c939SJed Brown ierr = DMCreateColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 181ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim); 18247c6ae99SBarry Smith if (isBAIJ) { 18347c6ae99SBarry Smith dd->w = nc; 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 PetscFunctionReturn(0); 19047c6ae99SBarry Smith } 19147c6ae99SBarry Smith 19247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 19347c6ae99SBarry Smith 194e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 19547c6ae99SBarry Smith { 19647c6ae99SBarry Smith PetscErrorCode ierr; 19747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col; 19847c6ae99SBarry Smith PetscInt ncolors; 19947c6ae99SBarry Smith MPI_Comm comm; 200bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by; 201aa219208SBarry Smith DMDAStencilType st; 20247c6ae99SBarry Smith ISColoringValue *colors; 20347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 20447c6ae99SBarry Smith 20547c6ae99SBarry Smith PetscFunctionBegin; 20647c6ae99SBarry Smith /* 20747c6ae99SBarry Smith nc - number of components per grid point 20847c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 20947c6ae99SBarry Smith 21047c6ae99SBarry Smith */ 2111321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 21247c6ae99SBarry Smith col = 2*s + 1; 213aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 214aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 21547c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 21647c6ae99SBarry Smith 21747c6ae99SBarry Smith /* special case as taught to us by Paul Hovland */ 218aa219208SBarry Smith if (st == DMDA_STENCIL_STAR && s == 1) { 219e727c939SJed Brown ierr = DMCreateColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr); 22047c6ae99SBarry Smith } else { 22147c6ae99SBarry Smith 222bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\ 22347c6ae99SBarry Smith by 2*stencil_width + 1 (%d)\n", m, col); 224bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC && (n % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\ 22547c6ae99SBarry Smith by 2*stencil_width + 1 (%d)\n", n, col); 22647c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 22747c6ae99SBarry Smith if (!dd->localcoloring) { 228785e854fSJed Brown ierr = PetscMalloc1(nc*nx*ny,&colors);CHKERRQ(ierr); 22947c6ae99SBarry Smith ii = 0; 23047c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 23147c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 23247c6ae99SBarry Smith for (k=0; k<nc; k++) { 23347c6ae99SBarry Smith colors[ii++] = k + nc*((i % col) + col*(j % col)); 23447c6ae99SBarry Smith } 23547c6ae99SBarry Smith } 23647c6ae99SBarry Smith } 23747c6ae99SBarry Smith ncolors = nc + nc*(col-1 + col*(col-1)); 238aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,PETSC_OWN_POINTER,&dd->localcoloring);CHKERRQ(ierr); 23947c6ae99SBarry Smith } 24047c6ae99SBarry Smith *coloring = dd->localcoloring; 2415bdb020cSBarry Smith } else if (ctype == IS_COLORING_LOCAL) { 24247c6ae99SBarry Smith if (!dd->ghostedcoloring) { 243785e854fSJed Brown ierr = PetscMalloc1(nc*gnx*gny,&colors);CHKERRQ(ierr); 24447c6ae99SBarry Smith ii = 0; 24547c6ae99SBarry Smith for (j=gys; j<gys+gny; j++) { 24647c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 24747c6ae99SBarry Smith for (k=0; k<nc; k++) { 24847c6ae99SBarry Smith /* the complicated stuff is to handle periodic boundaries */ 24947c6ae99SBarry Smith colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col)); 25047c6ae99SBarry Smith } 25147c6ae99SBarry Smith } 25247c6ae99SBarry Smith } 25347c6ae99SBarry Smith ncolors = nc + nc*(col - 1 + col*(col-1)); 254aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,PETSC_OWN_POINTER,&dd->ghostedcoloring);CHKERRQ(ierr); 25547c6ae99SBarry Smith /* PetscIntView(ncolors,(PetscInt*)colors,0); */ 25647c6ae99SBarry Smith 2575bdb020cSBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_LOCAL);CHKERRQ(ierr); 25847c6ae99SBarry Smith } 25947c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 260ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 26147c6ae99SBarry Smith } 26247c6ae99SBarry Smith ierr = ISColoringReference(*coloring);CHKERRQ(ierr); 26347c6ae99SBarry Smith PetscFunctionReturn(0); 26447c6ae99SBarry Smith } 26547c6ae99SBarry Smith 26647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 26747c6ae99SBarry Smith 268e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 26947c6ae99SBarry Smith { 27047c6ae99SBarry Smith PetscErrorCode ierr; 27147c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,gxs,gys,gnx,gny,m,n,p,dim,s,k,nc,col,zs,gzs,ii,l,nz,gnz,M,N,P; 27247c6ae99SBarry Smith PetscInt ncolors; 27347c6ae99SBarry Smith MPI_Comm comm; 274bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by,bz; 275aa219208SBarry Smith DMDAStencilType st; 27647c6ae99SBarry Smith ISColoringValue *colors; 27747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 27847c6ae99SBarry Smith 27947c6ae99SBarry Smith PetscFunctionBegin; 28047c6ae99SBarry Smith /* 28147c6ae99SBarry Smith nc - number of components per grid point 28247c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 28347c6ae99SBarry Smith 28447c6ae99SBarry Smith */ 2851321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 28647c6ae99SBarry Smith col = 2*s + 1; 287bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\ 28847c6ae99SBarry Smith by 2*stencil_width + 1\n"); 289bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\ 29047c6ae99SBarry Smith by 2*stencil_width + 1\n"); 291bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\ 29247c6ae99SBarry Smith by 2*stencil_width + 1\n"); 29347c6ae99SBarry Smith 294aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 295aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 29647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 29747c6ae99SBarry Smith 29847c6ae99SBarry Smith /* create the coloring */ 29947c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 30047c6ae99SBarry Smith if (!dd->localcoloring) { 301785e854fSJed Brown ierr = PetscMalloc1(nc*nx*ny*nz,&colors);CHKERRQ(ierr); 30247c6ae99SBarry Smith ii = 0; 30347c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 30447c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 30547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 30647c6ae99SBarry Smith for (l=0; l<nc; l++) { 30747c6ae99SBarry Smith colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col)); 30847c6ae99SBarry Smith } 30947c6ae99SBarry Smith } 31047c6ae99SBarry Smith } 31147c6ae99SBarry Smith } 31247c6ae99SBarry Smith ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1)); 313aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,PETSC_OWN_POINTER,&dd->localcoloring);CHKERRQ(ierr); 31447c6ae99SBarry Smith } 31547c6ae99SBarry Smith *coloring = dd->localcoloring; 3165bdb020cSBarry Smith } else if (ctype == IS_COLORING_LOCAL) { 31747c6ae99SBarry Smith if (!dd->ghostedcoloring) { 318785e854fSJed Brown ierr = PetscMalloc1(nc*gnx*gny*gnz,&colors);CHKERRQ(ierr); 31947c6ae99SBarry Smith ii = 0; 32047c6ae99SBarry Smith for (k=gzs; k<gzs+gnz; k++) { 32147c6ae99SBarry Smith for (j=gys; j<gys+gny; j++) { 32247c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 32347c6ae99SBarry Smith for (l=0; l<nc; l++) { 32447c6ae99SBarry Smith /* the complicated stuff is to handle periodic boundaries */ 32547c6ae99SBarry Smith colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col)); 32647c6ae99SBarry Smith } 32747c6ae99SBarry Smith } 32847c6ae99SBarry Smith } 32947c6ae99SBarry Smith } 33047c6ae99SBarry Smith ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1)); 331aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,PETSC_OWN_POINTER,&dd->ghostedcoloring);CHKERRQ(ierr); 3325bdb020cSBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_LOCAL);CHKERRQ(ierr); 33347c6ae99SBarry Smith } 33447c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 335ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 33647c6ae99SBarry Smith ierr = ISColoringReference(*coloring);CHKERRQ(ierr); 33747c6ae99SBarry Smith PetscFunctionReturn(0); 33847c6ae99SBarry Smith } 33947c6ae99SBarry Smith 34047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 34147c6ae99SBarry Smith 342e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 34347c6ae99SBarry Smith { 34447c6ae99SBarry Smith PetscErrorCode ierr; 34547c6ae99SBarry Smith PetscInt xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col; 34647c6ae99SBarry Smith PetscInt ncolors; 34747c6ae99SBarry Smith MPI_Comm comm; 348bff4a2f0SMatthew G. Knepley DMBoundaryType bx; 34947c6ae99SBarry Smith ISColoringValue *colors; 35047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 35147c6ae99SBarry Smith 35247c6ae99SBarry Smith PetscFunctionBegin; 35347c6ae99SBarry Smith /* 35447c6ae99SBarry Smith nc - number of components per grid point 35547c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 35647c6ae99SBarry Smith 35747c6ae99SBarry Smith */ 3581321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr); 35947c6ae99SBarry Smith col = 2*s + 1; 36047c6ae99SBarry Smith 361bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points %d is divisible\n\ 36231e6f798SBarry Smith by 2*stencil_width + 1 %d\n",(int)m,(int)col); 36347c6ae99SBarry Smith 364aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr); 365aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr); 36647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 36747c6ae99SBarry Smith 36847c6ae99SBarry Smith /* create the coloring */ 36947c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 37047c6ae99SBarry Smith if (!dd->localcoloring) { 371785e854fSJed Brown ierr = PetscMalloc1(nc*nx,&colors);CHKERRQ(ierr); 372ae4f298aSBarry Smith if (dd->ofillcols) { 373ae4f298aSBarry Smith PetscInt tc = 0; 374ae4f298aSBarry Smith for (i=0; i<nc; i++) tc += (PetscInt) (dd->ofillcols[i] > 0); 375ae4f298aSBarry Smith i1 = 0; 376ae4f298aSBarry Smith for (i=xs; i<xs+nx; i++) { 377ae4f298aSBarry Smith for (l=0; l<nc; l++) { 378ae4f298aSBarry Smith if (dd->ofillcols[l] && (i % col)) { 379ae4f298aSBarry Smith colors[i1++] = nc - 1 + tc*((i % col) - 1) + dd->ofillcols[l]; 380ae4f298aSBarry Smith } else { 381ae4f298aSBarry Smith colors[i1++] = l; 382ae4f298aSBarry Smith } 383ae4f298aSBarry Smith } 384ae4f298aSBarry Smith } 385ae4f298aSBarry Smith ncolors = nc + 2*s*tc; 386ae4f298aSBarry Smith } else { 38747c6ae99SBarry Smith i1 = 0; 38847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 38947c6ae99SBarry Smith for (l=0; l<nc; l++) { 39047c6ae99SBarry Smith colors[i1++] = l + nc*(i % col); 39147c6ae99SBarry Smith } 39247c6ae99SBarry Smith } 39347c6ae99SBarry Smith ncolors = nc + nc*(col-1); 394ae4f298aSBarry Smith } 395aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,PETSC_OWN_POINTER,&dd->localcoloring);CHKERRQ(ierr); 39647c6ae99SBarry Smith } 39747c6ae99SBarry Smith *coloring = dd->localcoloring; 3985bdb020cSBarry Smith } else if (ctype == IS_COLORING_LOCAL) { 39947c6ae99SBarry Smith if (!dd->ghostedcoloring) { 400785e854fSJed Brown ierr = PetscMalloc1(nc*gnx,&colors);CHKERRQ(ierr); 40147c6ae99SBarry Smith i1 = 0; 40247c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 40347c6ae99SBarry Smith for (l=0; l<nc; l++) { 40447c6ae99SBarry Smith /* the complicated stuff is to handle periodic boundaries */ 40547c6ae99SBarry Smith colors[i1++] = l + nc*(SetInRange(i,m) % col); 40647c6ae99SBarry Smith } 40747c6ae99SBarry Smith } 40847c6ae99SBarry Smith ncolors = nc + nc*(col-1); 409aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,PETSC_OWN_POINTER,&dd->ghostedcoloring);CHKERRQ(ierr); 4105bdb020cSBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_LOCAL);CHKERRQ(ierr); 41147c6ae99SBarry Smith } 41247c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 413ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 41447c6ae99SBarry Smith ierr = ISColoringReference(*coloring);CHKERRQ(ierr); 41547c6ae99SBarry Smith PetscFunctionReturn(0); 41647c6ae99SBarry Smith } 41747c6ae99SBarry Smith 418e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring) 41947c6ae99SBarry Smith { 42047c6ae99SBarry Smith PetscErrorCode ierr; 42147c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc; 42247c6ae99SBarry Smith PetscInt ncolors; 42347c6ae99SBarry Smith MPI_Comm comm; 424bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by; 42547c6ae99SBarry Smith ISColoringValue *colors; 42647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 42747c6ae99SBarry Smith 42847c6ae99SBarry Smith PetscFunctionBegin; 42947c6ae99SBarry Smith /* 43047c6ae99SBarry Smith nc - number of components per grid point 43147c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 43247c6ae99SBarry Smith 43347c6ae99SBarry Smith */ 4341321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr); 435aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 436aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 43747c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 43847c6ae99SBarry Smith 439bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && (m % 5)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible by 5\n"); 440bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC && (n % 5)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible by 5\n"); 44147c6ae99SBarry Smith 44247c6ae99SBarry Smith /* create the coloring */ 44347c6ae99SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 44447c6ae99SBarry Smith if (!dd->localcoloring) { 445785e854fSJed Brown ierr = PetscMalloc1(nc*nx*ny,&colors);CHKERRQ(ierr); 44647c6ae99SBarry Smith ii = 0; 44747c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 44847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 44947c6ae99SBarry Smith for (k=0; k<nc; k++) { 45047c6ae99SBarry Smith colors[ii++] = k + nc*((3*j+i) % 5); 45147c6ae99SBarry Smith } 45247c6ae99SBarry Smith } 45347c6ae99SBarry Smith } 45447c6ae99SBarry Smith ncolors = 5*nc; 455aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,PETSC_OWN_POINTER,&dd->localcoloring);CHKERRQ(ierr); 45647c6ae99SBarry Smith } 45747c6ae99SBarry Smith *coloring = dd->localcoloring; 4585bdb020cSBarry Smith } else if (ctype == IS_COLORING_LOCAL) { 45947c6ae99SBarry Smith if (!dd->ghostedcoloring) { 460785e854fSJed Brown ierr = PetscMalloc1(nc*gnx*gny,&colors);CHKERRQ(ierr); 46147c6ae99SBarry Smith ii = 0; 46247c6ae99SBarry Smith for (j=gys; j<gys+gny; j++) { 46347c6ae99SBarry Smith for (i=gxs; i<gxs+gnx; i++) { 46447c6ae99SBarry Smith for (k=0; k<nc; k++) { 46547c6ae99SBarry Smith colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5); 46647c6ae99SBarry Smith } 46747c6ae99SBarry Smith } 46847c6ae99SBarry Smith } 46947c6ae99SBarry Smith ncolors = 5*nc; 470aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,PETSC_OWN_POINTER,&dd->ghostedcoloring);CHKERRQ(ierr); 4715bdb020cSBarry Smith ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_LOCAL);CHKERRQ(ierr); 47247c6ae99SBarry Smith } 47347c6ae99SBarry Smith *coloring = dd->ghostedcoloring; 474ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype); 47547c6ae99SBarry Smith PetscFunctionReturn(0); 47647c6ae99SBarry Smith } 47747c6ae99SBarry Smith 47847c6ae99SBarry Smith /* =========================================================================== */ 479950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM,Mat); 480ce308e1dSBarry Smith extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM,Mat); 481950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM,Mat); 482950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM,Mat); 483950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM,Mat); 484950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM,Mat); 485950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM,Mat); 486950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM,Mat); 487950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM,Mat); 488950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM,Mat); 48947c6ae99SBarry Smith 4908bbdbebaSMatthew G Knepley /*@C 491c688c046SMatthew G Knepley MatSetupDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix 49247c6ae99SBarry Smith 49347c6ae99SBarry Smith Logically Collective on Mat 49447c6ae99SBarry Smith 49547c6ae99SBarry Smith Input Parameters: 49647c6ae99SBarry Smith + mat - the matrix 49747c6ae99SBarry Smith - da - the da 49847c6ae99SBarry Smith 49947c6ae99SBarry Smith Level: intermediate 50047c6ae99SBarry Smith 50147c6ae99SBarry Smith @*/ 502c688c046SMatthew G Knepley PetscErrorCode MatSetupDM(Mat mat,DM da) 50347c6ae99SBarry Smith { 50447c6ae99SBarry Smith PetscErrorCode ierr; 50547c6ae99SBarry Smith 50647c6ae99SBarry Smith PetscFunctionBegin; 50747c6ae99SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 50847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 509c688c046SMatthew G Knepley ierr = PetscTryMethod(mat,"MatSetupDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr); 51047c6ae99SBarry Smith PetscFunctionReturn(0); 51147c6ae99SBarry Smith } 51247c6ae99SBarry Smith 5137087cfbeSBarry Smith PetscErrorCode MatView_MPI_DA(Mat A,PetscViewer viewer) 51447c6ae99SBarry Smith { 5159a42bb27SBarry Smith DM da; 51647c6ae99SBarry Smith PetscErrorCode ierr; 51747c6ae99SBarry Smith const char *prefix; 51847c6ae99SBarry Smith Mat Anatural; 51947c6ae99SBarry Smith AO ao; 52047c6ae99SBarry Smith PetscInt rstart,rend,*petsc,i; 52147c6ae99SBarry Smith IS is; 52247c6ae99SBarry Smith MPI_Comm comm; 52374388724SJed Brown PetscViewerFormat format; 52447c6ae99SBarry Smith 52547c6ae99SBarry Smith PetscFunctionBegin; 52674388724SJed Brown /* Check whether we are just printing info, in which case MatView() already viewed everything we wanted to view */ 52774388724SJed Brown ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 52874388724SJed Brown if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0); 52974388724SJed Brown 53047c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 531c688c046SMatthew G Knepley ierr = MatGetDM(A, &da);CHKERRQ(ierr); 532ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA"); 53347c6ae99SBarry Smith 534aa219208SBarry Smith ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr); 53547c6ae99SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 536854ce69bSBarry Smith ierr = PetscMalloc1(rend-rstart,&petsc);CHKERRQ(ierr); 53747c6ae99SBarry Smith for (i=rstart; i<rend; i++) petsc[i-rstart] = i; 53847c6ae99SBarry Smith ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr); 53947c6ae99SBarry Smith ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 54047c6ae99SBarry Smith 54147c6ae99SBarry Smith /* call viewer on natural ordering */ 5427dae84e0SHong Zhang ierr = MatCreateSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr); 543fcfd50ebSBarry Smith ierr = ISDestroy(&is);CHKERRQ(ierr); 54447c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr); 54547c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr); 54647c6ae99SBarry Smith ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr); 547539c167fSBarry Smith ierr = MatView(Anatural,viewer);CHKERRQ(ierr); 548fcfd50ebSBarry Smith ierr = MatDestroy(&Anatural);CHKERRQ(ierr); 54947c6ae99SBarry Smith PetscFunctionReturn(0); 55047c6ae99SBarry Smith } 55147c6ae99SBarry Smith 5527087cfbeSBarry Smith PetscErrorCode MatLoad_MPI_DA(Mat A,PetscViewer viewer) 55347c6ae99SBarry Smith { 5549a42bb27SBarry Smith DM da; 55547c6ae99SBarry Smith PetscErrorCode ierr; 55647c6ae99SBarry Smith Mat Anatural,Aapp; 55747c6ae99SBarry Smith AO ao; 558539c167fSBarry Smith PetscInt rstart,rend,*app,i,m,n,M,N; 55947c6ae99SBarry Smith IS is; 56047c6ae99SBarry Smith MPI_Comm comm; 56147c6ae99SBarry Smith 56247c6ae99SBarry Smith PetscFunctionBegin; 56347c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 564c688c046SMatthew G Knepley ierr = MatGetDM(A, &da);CHKERRQ(ierr); 565ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA"); 56647c6ae99SBarry Smith 56747c6ae99SBarry Smith /* Load the matrix in natural ordering */ 568ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&Anatural);CHKERRQ(ierr); 56947c6ae99SBarry Smith ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr); 570539c167fSBarry Smith ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 571539c167fSBarry Smith ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 572539c167fSBarry Smith ierr = MatSetSizes(Anatural,m,n,M,N);CHKERRQ(ierr); 57347c6ae99SBarry Smith ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr); 57447c6ae99SBarry Smith 57547c6ae99SBarry Smith /* Map natural ordering to application ordering and create IS */ 576aa219208SBarry Smith ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr); 57747c6ae99SBarry Smith ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr); 578854ce69bSBarry Smith ierr = PetscMalloc1(rend-rstart,&app);CHKERRQ(ierr); 57947c6ae99SBarry Smith for (i=rstart; i<rend; i++) app[i-rstart] = i; 58047c6ae99SBarry Smith ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr); 58147c6ae99SBarry Smith ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 58247c6ae99SBarry Smith 58347c6ae99SBarry Smith /* Do permutation and replace header */ 5847dae84e0SHong Zhang ierr = MatCreateSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr); 58528be2f97SBarry Smith ierr = MatHeaderReplace(A,&Aapp);CHKERRQ(ierr); 586fcfd50ebSBarry Smith ierr = ISDestroy(&is);CHKERRQ(ierr); 587fcfd50ebSBarry Smith ierr = MatDestroy(&Anatural);CHKERRQ(ierr); 58847c6ae99SBarry Smith PetscFunctionReturn(0); 58947c6ae99SBarry Smith } 59047c6ae99SBarry Smith 591b412c318SBarry Smith PetscErrorCode DMCreateMatrix_DA(DM da, Mat *J) 59247c6ae99SBarry Smith { 59347c6ae99SBarry Smith PetscErrorCode ierr; 59447c6ae99SBarry Smith PetscInt dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P; 59547c6ae99SBarry Smith Mat A; 59647c6ae99SBarry Smith MPI_Comm comm; 59719fd82e9SBarry Smith MatType Atype; 59837d0c07bSMatthew G Knepley PetscSection section, sectionGlobal; 5990298fd71SBarry Smith void (*aij)(void)=NULL,(*baij)(void)=NULL,(*sbaij)(void)=NULL; 600b412c318SBarry Smith MatType mtype; 60147c6ae99SBarry Smith PetscMPIInt size; 60247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 60347c6ae99SBarry Smith 60447c6ae99SBarry Smith PetscFunctionBegin; 605607a6623SBarry Smith ierr = MatInitializePackage();CHKERRQ(ierr); 606b412c318SBarry Smith mtype = da->mattype; 60747c6ae99SBarry Smith 60837d0c07bSMatthew G Knepley ierr = DMGetDefaultSection(da, §ion);CHKERRQ(ierr); 60937d0c07bSMatthew G Knepley if (section) { 61037d0c07bSMatthew G Knepley PetscInt bs = -1; 61137d0c07bSMatthew G Knepley PetscInt localSize; 61237d0c07bSMatthew G Knepley PetscBool isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isSymmetric; 61337d0c07bSMatthew G Knepley 61437d0c07bSMatthew G Knepley ierr = DMGetDefaultGlobalSection(da, §ionGlobal);CHKERRQ(ierr); 61537d0c07bSMatthew G Knepley ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);CHKERRQ(ierr); 61682f516ccSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)da), J);CHKERRQ(ierr); 61737d0c07bSMatthew G Knepley ierr = MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 61837d0c07bSMatthew G Knepley ierr = MatSetType(*J, mtype);CHKERRQ(ierr); 61937d0c07bSMatthew G Knepley ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 62037d0c07bSMatthew G Knepley ierr = PetscStrcmp(mtype, MATSHELL, &isShell);CHKERRQ(ierr); 62137d0c07bSMatthew G Knepley ierr = PetscStrcmp(mtype, MATBAIJ, &isBlock);CHKERRQ(ierr); 62237d0c07bSMatthew G Knepley ierr = PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock);CHKERRQ(ierr); 62337d0c07bSMatthew G Knepley ierr = PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock);CHKERRQ(ierr); 62437d0c07bSMatthew G Knepley ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr); 62537d0c07bSMatthew G Knepley ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr); 62637d0c07bSMatthew G Knepley ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr); 62737d0c07bSMatthew G Knepley /* Check for symmetric storage */ 62837d0c07bSMatthew G Knepley isSymmetric = (PetscBool) (isSymBlock || isSymSeqBlock || isSymMPIBlock); 62937d0c07bSMatthew G Knepley if (isSymmetric) { 63037d0c07bSMatthew G Knepley ierr = MatSetOption(*J, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr); 63137d0c07bSMatthew G Knepley } 63237d0c07bSMatthew G Knepley if (!isShell) { 63337d0c07bSMatthew G Knepley PetscInt *dnz, *onz, *dnzu, *onzu, bsLocal; 63437d0c07bSMatthew G Knepley 63537d0c07bSMatthew G Knepley if (bs < 0) { 63637d0c07bSMatthew G Knepley if (isBlock || isSeqBlock || isMPIBlock || isSymBlock || isSymSeqBlock || isSymMPIBlock) { 63737d0c07bSMatthew G Knepley PetscInt pStart, pEnd, p, dof; 63837d0c07bSMatthew G Knepley 63937d0c07bSMatthew G Knepley ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr); 64037d0c07bSMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 64137d0c07bSMatthew G Knepley ierr = PetscSectionGetDof(sectionGlobal, p, &dof);CHKERRQ(ierr); 64237d0c07bSMatthew G Knepley if (dof) { 64337d0c07bSMatthew G Knepley bs = dof; 64437d0c07bSMatthew G Knepley break; 64537d0c07bSMatthew G Knepley } 64637d0c07bSMatthew G Knepley } 64737d0c07bSMatthew G Knepley } else { 64837d0c07bSMatthew G Knepley bs = 1; 64937d0c07bSMatthew G Knepley } 65037d0c07bSMatthew G Knepley /* Must have same blocksize on all procs (some might have no points) */ 65137d0c07bSMatthew G Knepley bsLocal = bs; 652b2566f29SBarry Smith ierr = MPIU_Allreduce(&bsLocal, &bs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 65337d0c07bSMatthew G Knepley } 6541795a4d1SJed Brown ierr = PetscCalloc4(localSize/bs, &dnz, localSize/bs, &onz, localSize/bs, &dnzu, localSize/bs, &onzu);CHKERRQ(ierr); 655552f7358SJed Brown /* ierr = DMPlexPreallocateOperator(dm, bs, section, sectionGlobal, dnz, onz, dnzu, onzu, *J, fillMatrix);CHKERRQ(ierr); */ 65637d0c07bSMatthew G Knepley ierr = PetscFree4(dnz, onz, dnzu, onzu);CHKERRQ(ierr); 65737d0c07bSMatthew G Knepley } 65837d0c07bSMatthew G Knepley } 65947c6ae99SBarry Smith /* 66047c6ae99SBarry Smith m 66147c6ae99SBarry Smith ------------------------------------------------------ 66247c6ae99SBarry Smith | | 66347c6ae99SBarry Smith | | 66447c6ae99SBarry Smith | ---------------------- | 66547c6ae99SBarry Smith | | | | 66647c6ae99SBarry Smith n | ny | | | 66747c6ae99SBarry Smith | | | | 66847c6ae99SBarry Smith | .--------------------- | 66947c6ae99SBarry Smith | (xs,ys) nx | 67047c6ae99SBarry Smith | . | 67147c6ae99SBarry Smith | (gxs,gys) | 67247c6ae99SBarry Smith | | 67347c6ae99SBarry Smith ----------------------------------------------------- 67447c6ae99SBarry Smith */ 67547c6ae99SBarry Smith 67647c6ae99SBarry Smith /* 67747c6ae99SBarry Smith nc - number of components per grid point 67847c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 67947c6ae99SBarry Smith 68047c6ae99SBarry Smith */ 681e30e807fSPeter Brune M = dd->M; 682e30e807fSPeter Brune N = dd->N; 683e30e807fSPeter Brune P = dd->P; 684c73cfb54SMatthew G. Knepley dim = da->dim; 685e30e807fSPeter Brune dof = dd->w; 686e30e807fSPeter Brune /* ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); */ 687aa219208SBarry Smith ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr); 68847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 68947c6ae99SBarry Smith ierr = MatCreate(comm,&A);CHKERRQ(ierr); 69047c6ae99SBarry Smith ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr); 691b412c318SBarry Smith ierr = MatSetType(A,mtype);CHKERRQ(ierr); 69295ee5b0eSBarry Smith ierr = MatSetDM(A,da);CHKERRQ(ierr); 693*b06ff27eSHong Zhang if (da->structure_only) { 694*b06ff27eSHong Zhang ierr = MatSetOption(A,MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr); 695*b06ff27eSHong Zhang } 69647c6ae99SBarry Smith ierr = MatSetFromOptions(A);CHKERRQ(ierr); 69747c6ae99SBarry Smith ierr = MatGetType(A,&Atype);CHKERRQ(ierr); 69847c6ae99SBarry Smith /* 699aa219208SBarry Smith We do not provide a getmatrix function in the DMDA operations because 700aa219208SBarry Smith the basic DMDA does not know about matrices. We think of DMDA as being more 70147c6ae99SBarry Smith more low-level than matrices. This is kind of cheating but, cause sometimes 702aa219208SBarry Smith we think of DMDA has higher level than matrices. 70347c6ae99SBarry Smith 70447c6ae99SBarry Smith We could switch based on Atype (or mtype), but we do not since the 70547c6ae99SBarry Smith specialized setting routines depend only the particular preallocation 70647c6ae99SBarry Smith details of the matrix, not the type itself. 70747c6ae99SBarry Smith */ 70847c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 70947c6ae99SBarry Smith if (!aij) { 71047c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 71147c6ae99SBarry Smith } 71247c6ae99SBarry Smith if (!aij) { 71347c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr); 71447c6ae99SBarry Smith if (!baij) { 71547c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr); 71647c6ae99SBarry Smith } 71747c6ae99SBarry Smith if (!baij) { 71847c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr); 71947c6ae99SBarry Smith if (!sbaij) { 72047c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr); 72147c6ae99SBarry Smith } 72247c6ae99SBarry Smith } 72347c6ae99SBarry Smith } 72447c6ae99SBarry Smith if (aij) { 72547c6ae99SBarry Smith if (dim == 1) { 726ce308e1dSBarry Smith if (dd->ofill) { 727ce308e1dSBarry Smith ierr = DMCreateMatrix_DA_1d_MPIAIJ_Fill(da,A);CHKERRQ(ierr); 728ce308e1dSBarry Smith } else { 729950540a4SJed Brown ierr = DMCreateMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr); 730ce308e1dSBarry Smith } 73147c6ae99SBarry Smith } else if (dim == 2) { 73247c6ae99SBarry Smith if (dd->ofill) { 733950540a4SJed Brown ierr = DMCreateMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr); 73447c6ae99SBarry Smith } else { 735950540a4SJed Brown ierr = DMCreateMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr); 73647c6ae99SBarry Smith } 73747c6ae99SBarry Smith } else if (dim == 3) { 73847c6ae99SBarry Smith if (dd->ofill) { 739950540a4SJed Brown ierr = DMCreateMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr); 74047c6ae99SBarry Smith } else { 741950540a4SJed Brown ierr = DMCreateMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr); 74247c6ae99SBarry Smith } 74347c6ae99SBarry Smith } 74447c6ae99SBarry Smith } else if (baij) { 74547c6ae99SBarry Smith if (dim == 2) { 746950540a4SJed Brown ierr = DMCreateMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr); 74747c6ae99SBarry Smith } else if (dim == 3) { 748950540a4SJed Brown ierr = DMCreateMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr); 749ce94432eSBarry Smith } else SETERRQ3(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension! Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim); 75047c6ae99SBarry Smith } else if (sbaij) { 75147c6ae99SBarry Smith if (dim == 2) { 752950540a4SJed Brown ierr = DMCreateMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr); 75347c6ae99SBarry Smith } else if (dim == 3) { 754950540a4SJed Brown ierr = DMCreateMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr); 755ce94432eSBarry Smith } else SETERRQ3(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension! Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim); 756869776cdSLisandro Dalcin } else { 75745b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 758869776cdSLisandro Dalcin ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 7592949035bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 760869776cdSLisandro Dalcin ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr); 76147c6ae99SBarry Smith } 762aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr); 76347c6ae99SBarry Smith ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr); 764c688c046SMatthew G Knepley ierr = MatSetDM(A,da);CHKERRQ(ierr); 76547c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 76647c6ae99SBarry Smith if (size > 1) { 76747c6ae99SBarry Smith /* change viewer to display matrix in natural ordering */ 76847c6ae99SBarry Smith ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MPI_DA);CHKERRQ(ierr); 76947c6ae99SBarry Smith ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void))MatLoad_MPI_DA);CHKERRQ(ierr); 77047c6ae99SBarry Smith } 77147c6ae99SBarry Smith *J = A; 77247c6ae99SBarry Smith PetscFunctionReturn(0); 77347c6ae99SBarry Smith } 77447c6ae99SBarry Smith 77547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 776950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM da,Mat J) 77747c6ae99SBarry Smith { 77847c6ae99SBarry Smith PetscErrorCode ierr; 779c1154cd5SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny,m,n,dim,s,*cols = NULL,k,nc,*rows = NULL,col,cnt,l,p,M,N; 78047c6ae99SBarry Smith PetscInt lstart,lend,pstart,pend,*dnz,*onz; 78147c6ae99SBarry Smith MPI_Comm comm; 78247c6ae99SBarry Smith PetscScalar *values; 783bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by; 78445b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 785aa219208SBarry Smith DMDAStencilType st; 786c1154cd5SBarry Smith PetscBool removedups = PETSC_FALSE; 78747c6ae99SBarry Smith 78847c6ae99SBarry Smith PetscFunctionBegin; 78947c6ae99SBarry Smith /* 79047c6ae99SBarry Smith nc - number of components per grid point 79147c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 79247c6ae99SBarry Smith 79347c6ae99SBarry Smith */ 794c1154cd5SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,&M,&N,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 79547c6ae99SBarry Smith col = 2*s + 1; 796c1154cd5SBarry Smith /* 797c1154cd5SBarry Smith With one processor in periodic domains in a skinny dimension the code will label nonzero columns multiple times 798c1154cd5SBarry Smith because of "wrapping" around the end of the domain hitting an entry already counted in the other direction. 799c1154cd5SBarry Smith */ 800c1154cd5SBarry Smith if (M == 1 && 2*s >= m) removedups = PETSC_TRUE; 801c1154cd5SBarry Smith if (N == 1 && 2*s >= n) removedups = PETSC_TRUE; 802aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 803aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 80447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 80547c6ae99SBarry Smith 806dcca6d9dSJed Brown ierr = PetscMalloc2(nc,&rows,col*col*nc*nc,&cols);CHKERRQ(ierr); 8071411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 80847c6ae99SBarry Smith 80906ca8cadSBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 81047c6ae99SBarry Smith /* determine the matrix preallocation information */ 81147c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr); 81247c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 81347c6ae99SBarry Smith 814bff4a2f0SMatthew G. Knepley pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 815bff4a2f0SMatthew G. Knepley pend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 81647c6ae99SBarry Smith 81747c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 81847c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 81947c6ae99SBarry Smith 820bff4a2f0SMatthew G. Knepley lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 821bff4a2f0SMatthew G. Knepley lend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 82247c6ae99SBarry Smith 82347c6ae99SBarry Smith cnt = 0; 82447c6ae99SBarry Smith for (k=0; k<nc; k++) { 82547c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 82647c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 827aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */ 82847c6ae99SBarry Smith cols[cnt++] = k + nc*(slot + gnx*l + p); 82947c6ae99SBarry Smith } 83047c6ae99SBarry Smith } 83147c6ae99SBarry Smith } 83247c6ae99SBarry Smith rows[k] = k + nc*(slot); 83347c6ae99SBarry Smith } 834c1154cd5SBarry Smith if (removedups) { 835c1154cd5SBarry Smith ierr = MatPreallocateSetLocalRemoveDups(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 836c1154cd5SBarry Smith } else { 837784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 83847c6ae99SBarry Smith } 83947c6ae99SBarry Smith } 840c1154cd5SBarry Smith } 841f73d5cc4SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 84247c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 84347c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 84447c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 84547c6ae99SBarry Smith 846784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 84747c6ae99SBarry Smith 84847c6ae99SBarry Smith /* 84947c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 85047c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 85147c6ae99SBarry Smith PETSc ordering. 85247c6ae99SBarry Smith */ 853fcfd50ebSBarry Smith if (!da->prealloc_only) { 8541795a4d1SJed Brown ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr); 85547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 85647c6ae99SBarry Smith 857bff4a2f0SMatthew G. Knepley pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 858bff4a2f0SMatthew G. Knepley pend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 85947c6ae99SBarry Smith 86047c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 86147c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 86247c6ae99SBarry Smith 863bff4a2f0SMatthew G. Knepley lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 864bff4a2f0SMatthew G. Knepley lend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 86547c6ae99SBarry Smith 86647c6ae99SBarry Smith cnt = 0; 86747c6ae99SBarry Smith for (k=0; k<nc; k++) { 86847c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 86947c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 870aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */ 87147c6ae99SBarry Smith cols[cnt++] = k + nc*(slot + gnx*l + p); 87247c6ae99SBarry Smith } 87347c6ae99SBarry Smith } 87447c6ae99SBarry Smith } 87547c6ae99SBarry Smith rows[k] = k + nc*(slot); 87647c6ae99SBarry Smith } 87747c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 87847c6ae99SBarry Smith } 87947c6ae99SBarry Smith } 88047c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 88147c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 88247c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 883189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 88447c6ae99SBarry Smith } 88547c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 88647c6ae99SBarry Smith PetscFunctionReturn(0); 88747c6ae99SBarry Smith } 88847c6ae99SBarry Smith 889950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J) 89047c6ae99SBarry Smith { 89147c6ae99SBarry Smith PetscErrorCode ierr; 89247c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 893c1154cd5SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,row,col,cnt,maxcnt = 0,l,p,M,N; 89447c6ae99SBarry Smith PetscInt lstart,lend,pstart,pend,*dnz,*onz; 89547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 89647c6ae99SBarry Smith PetscInt ifill_col,*ofill = dd->ofill, *dfill = dd->dfill; 89747c6ae99SBarry Smith MPI_Comm comm; 89847c6ae99SBarry Smith PetscScalar *values; 899bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by; 90045b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 901aa219208SBarry Smith DMDAStencilType st; 902c1154cd5SBarry Smith PetscBool removedups = PETSC_FALSE; 90347c6ae99SBarry Smith 90447c6ae99SBarry Smith PetscFunctionBegin; 90547c6ae99SBarry Smith /* 90647c6ae99SBarry Smith nc - number of components per grid point 90747c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 90847c6ae99SBarry Smith 90947c6ae99SBarry Smith */ 910c1154cd5SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,&M,&N,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 91147c6ae99SBarry Smith col = 2*s + 1; 912c1154cd5SBarry Smith /* 913c1154cd5SBarry Smith With one processor in periodic domains in a skinny dimension the code will label nonzero columns multiple times 914c1154cd5SBarry Smith because of "wrapping" around the end of the domain hitting an entry already counted in the other direction. 915c1154cd5SBarry Smith */ 916c1154cd5SBarry Smith if (M == 1 && 2*s >= m) removedups = PETSC_TRUE; 917c1154cd5SBarry Smith if (N == 1 && 2*s >= n) removedups = PETSC_TRUE; 918aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 919aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 92047c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 92147c6ae99SBarry Smith 9224b26d1cfSBarry Smith ierr = PetscMalloc1(col*col*nc,&cols);CHKERRQ(ierr); 9231411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 92447c6ae99SBarry Smith 92506ca8cadSBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 92647c6ae99SBarry Smith /* determine the matrix preallocation information */ 92747c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr); 92847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 92947c6ae99SBarry Smith 930bff4a2f0SMatthew G. Knepley pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 931bff4a2f0SMatthew G. Knepley pend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 93247c6ae99SBarry Smith 93347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 93447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 93547c6ae99SBarry Smith 936bff4a2f0SMatthew G. Knepley lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 937bff4a2f0SMatthew G. Knepley lend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 93847c6ae99SBarry Smith 93947c6ae99SBarry Smith for (k=0; k<nc; k++) { 94047c6ae99SBarry Smith cnt = 0; 94147c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 94247c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 94347c6ae99SBarry Smith if (l || p) { 944aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star */ 9458865f1eaSKarl Rupp for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p); 94647c6ae99SBarry Smith } 94747c6ae99SBarry Smith } else { 94847c6ae99SBarry Smith if (dfill) { 9498865f1eaSKarl Rupp for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p); 95047c6ae99SBarry Smith } else { 9518865f1eaSKarl Rupp for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p); 95247c6ae99SBarry Smith } 95347c6ae99SBarry Smith } 95447c6ae99SBarry Smith } 95547c6ae99SBarry Smith } 95647c6ae99SBarry Smith row = k + nc*(slot); 957c0ab637bSBarry Smith maxcnt = PetscMax(maxcnt,cnt); 958c1154cd5SBarry Smith if (removedups) { 959c1154cd5SBarry Smith ierr = MatPreallocateSetLocalRemoveDups(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 960c1154cd5SBarry Smith } else { 961784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 96247c6ae99SBarry Smith } 96347c6ae99SBarry Smith } 96447c6ae99SBarry Smith } 965c1154cd5SBarry Smith } 96647c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 96747c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 96847c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 969784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 97047c6ae99SBarry Smith 97147c6ae99SBarry Smith /* 97247c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 97347c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 97447c6ae99SBarry Smith PETSc ordering. 97547c6ae99SBarry Smith */ 976fcfd50ebSBarry Smith if (!da->prealloc_only) { 977c0ab637bSBarry Smith ierr = PetscCalloc1(maxcnt,&values);CHKERRQ(ierr); 97847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 97947c6ae99SBarry Smith 980bff4a2f0SMatthew G. Knepley pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 981bff4a2f0SMatthew G. Knepley pend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 98247c6ae99SBarry Smith 98347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 98447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 98547c6ae99SBarry Smith 986bff4a2f0SMatthew G. Knepley lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 987bff4a2f0SMatthew G. Knepley lend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 98847c6ae99SBarry Smith 98947c6ae99SBarry Smith for (k=0; k<nc; k++) { 99047c6ae99SBarry Smith cnt = 0; 99147c6ae99SBarry Smith for (l=lstart; l<lend+1; l++) { 99247c6ae99SBarry Smith for (p=pstart; p<pend+1; p++) { 99347c6ae99SBarry Smith if (l || p) { 994aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!l || !p)) { /* entries on star */ 9958865f1eaSKarl Rupp for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p); 99647c6ae99SBarry Smith } 99747c6ae99SBarry Smith } else { 99847c6ae99SBarry Smith if (dfill) { 9998865f1eaSKarl Rupp for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p); 100047c6ae99SBarry Smith } else { 10018865f1eaSKarl Rupp for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p); 100247c6ae99SBarry Smith } 100347c6ae99SBarry Smith } 100447c6ae99SBarry Smith } 100547c6ae99SBarry Smith } 100647c6ae99SBarry Smith row = k + nc*(slot); 100747c6ae99SBarry Smith ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 100847c6ae99SBarry Smith } 100947c6ae99SBarry Smith } 101047c6ae99SBarry Smith } 101147c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 101247c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 101347c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1014189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 101547c6ae99SBarry Smith } 101647c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 101747c6ae99SBarry Smith PetscFunctionReturn(0); 101847c6ae99SBarry Smith } 101947c6ae99SBarry Smith 102047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 102147c6ae99SBarry Smith 1022950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM da,Mat J) 102347c6ae99SBarry Smith { 102447c6ae99SBarry Smith PetscErrorCode ierr; 102547c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 10260298fd71SBarry Smith PetscInt m,n,dim,s,*cols = NULL,k,nc,*rows = NULL,col,cnt,l,p,*dnz = NULL,*onz = NULL; 1027c1154cd5SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk,M,N,P; 102847c6ae99SBarry Smith MPI_Comm comm; 102947c6ae99SBarry Smith PetscScalar *values; 1030bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by,bz; 103145b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 1032aa219208SBarry Smith DMDAStencilType st; 1033c1154cd5SBarry Smith PetscBool removedups = PETSC_FALSE; 103447c6ae99SBarry Smith 103547c6ae99SBarry Smith PetscFunctionBegin; 103647c6ae99SBarry Smith /* 103747c6ae99SBarry Smith nc - number of components per grid point 103847c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 103947c6ae99SBarry Smith 104047c6ae99SBarry Smith */ 1041c1154cd5SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 104247c6ae99SBarry Smith col = 2*s + 1; 104347c6ae99SBarry Smith 1044c1154cd5SBarry Smith /* 1045c1154cd5SBarry Smith With one processor in periodic domains in a skinny dimension the code will label nonzero columns multiple times 1046c1154cd5SBarry Smith because of "wrapping" around the end of the domain hitting an entry already counted in the other direction. 1047c1154cd5SBarry Smith */ 1048c1154cd5SBarry Smith if (M == 1 && 2*s >= m) removedups = PETSC_TRUE; 1049c1154cd5SBarry Smith if (N == 1 && 2*s >= n) removedups = PETSC_TRUE; 1050c1154cd5SBarry Smith if (P == 1 && 2*s >= p) removedups = PETSC_TRUE; 1051c1154cd5SBarry Smith 1052aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1053aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 105447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 105547c6ae99SBarry Smith 1056dcca6d9dSJed Brown ierr = PetscMalloc2(nc,&rows,col*col*col*nc*nc,&cols);CHKERRQ(ierr); 10571411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 105847c6ae99SBarry Smith 105906ca8cadSBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 106047c6ae99SBarry Smith /* determine the matrix preallocation information */ 106147c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr); 106247c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1063bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1064bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 106547c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1066bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1067bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 106847c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1069bff4a2f0SMatthew G. Knepley kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 1070bff4a2f0SMatthew G. Knepley kend = (bz == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 107147c6ae99SBarry Smith 107247c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 107347c6ae99SBarry Smith 107447c6ae99SBarry Smith cnt = 0; 107547c6ae99SBarry Smith for (l=0; l<nc; l++) { 107647c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 107747c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 107847c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1079aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 108047c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk); 108147c6ae99SBarry Smith } 108247c6ae99SBarry Smith } 108347c6ae99SBarry Smith } 108447c6ae99SBarry Smith } 108547c6ae99SBarry Smith rows[l] = l + nc*(slot); 108647c6ae99SBarry Smith } 1087c1154cd5SBarry Smith if (removedups) { 1088c1154cd5SBarry Smith ierr = MatPreallocateSetLocalRemoveDups(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 1089c1154cd5SBarry Smith } else { 1090784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 109147c6ae99SBarry Smith } 109247c6ae99SBarry Smith } 109347c6ae99SBarry Smith } 1094c1154cd5SBarry Smith } 1095f73d5cc4SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 109647c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 109747c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 109847c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1099784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 110047c6ae99SBarry Smith 110147c6ae99SBarry Smith /* 110247c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 110347c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 110447c6ae99SBarry Smith PETSc ordering. 110547c6ae99SBarry Smith */ 1106fcfd50ebSBarry Smith if (!da->prealloc_only) { 11071795a4d1SJed Brown ierr = PetscCalloc1(col*col*col*nc*nc*nc,&values);CHKERRQ(ierr); 110847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1109bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1110bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 111147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1112bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1113bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 111447c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1115bff4a2f0SMatthew G. Knepley kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 1116bff4a2f0SMatthew G. Knepley kend = (bz == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 111747c6ae99SBarry Smith 111847c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 111947c6ae99SBarry Smith 112047c6ae99SBarry Smith cnt = 0; 112147c6ae99SBarry Smith for (l=0; l<nc; l++) { 112247c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 112347c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 112447c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1125aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 112647c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk); 112747c6ae99SBarry Smith } 112847c6ae99SBarry Smith } 112947c6ae99SBarry Smith } 113047c6ae99SBarry Smith } 113147c6ae99SBarry Smith rows[l] = l + nc*(slot); 113247c6ae99SBarry Smith } 113347c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 113447c6ae99SBarry Smith } 113547c6ae99SBarry Smith } 113647c6ae99SBarry Smith } 113747c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 113847c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 113947c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1140189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 114147c6ae99SBarry Smith } 114247c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 114347c6ae99SBarry Smith PetscFunctionReturn(0); 114447c6ae99SBarry Smith } 114547c6ae99SBarry Smith 114647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 114747c6ae99SBarry Smith 1148ce308e1dSBarry Smith PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM da,Mat J) 1149ce308e1dSBarry Smith { 1150ce308e1dSBarry Smith PetscErrorCode ierr; 1151ce308e1dSBarry Smith DM_DA *dd = (DM_DA*)da->data; 1152ce308e1dSBarry Smith PetscInt xs,nx,i,j,gxs,gnx,row,k,l; 11538d4c968fSBarry Smith PetscInt m,dim,s,*cols = NULL,nc,cnt,maxcnt = 0,*ocols; 11540acb5bebSBarry Smith PetscInt *ofill = dd->ofill,*dfill = dd->dfill; 1155ce308e1dSBarry Smith PetscScalar *values; 1156bff4a2f0SMatthew G. Knepley DMBoundaryType bx; 115745b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 1158ce308e1dSBarry Smith PetscMPIInt rank,size; 1159ce308e1dSBarry Smith 1160ce308e1dSBarry Smith PetscFunctionBegin; 1161bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"With fill provided not implemented with periodic boundary conditions"); 1162ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr); 1163ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr); 1164ce308e1dSBarry Smith 1165ce308e1dSBarry Smith /* 1166ce308e1dSBarry Smith nc - number of components per grid point 1167ce308e1dSBarry Smith 1168ce308e1dSBarry Smith */ 1169ce308e1dSBarry Smith ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr); 1170ce308e1dSBarry Smith ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr); 1171ce308e1dSBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr); 1172ce308e1dSBarry Smith 1173ce308e1dSBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 11741795a4d1SJed Brown ierr = PetscCalloc2(nx*nc,&cols,nx*nc,&ocols);CHKERRQ(ierr); 1175ce308e1dSBarry Smith 1176ce308e1dSBarry Smith /* 1177ce308e1dSBarry Smith note should be smaller for first and last process with no periodic 1178ce308e1dSBarry Smith does not handle dfill 1179ce308e1dSBarry Smith */ 1180ce308e1dSBarry Smith cnt = 0; 1181ce308e1dSBarry Smith /* coupling with process to the left */ 1182ce308e1dSBarry Smith for (i=0; i<s; i++) { 1183ce308e1dSBarry Smith for (j=0; j<nc; j++) { 1184ce308e1dSBarry Smith ocols[cnt] = ((!rank) ? 0 : (s - i)*(ofill[j+1] - ofill[j])); 11850acb5bebSBarry Smith cols[cnt] = dfill[j+1] - dfill[j] + (s + i)*(ofill[j+1] - ofill[j]); 1186c0ab637bSBarry Smith maxcnt = PetscMax(maxcnt,ocols[cnt]+cols[cnt]); 1187ce308e1dSBarry Smith cnt++; 1188ce308e1dSBarry Smith } 1189ce308e1dSBarry Smith } 1190ce308e1dSBarry Smith for (i=s; i<nx-s; i++) { 1191ce308e1dSBarry Smith for (j=0; j<nc; j++) { 11920acb5bebSBarry Smith cols[cnt] = dfill[j+1] - dfill[j] + 2*s*(ofill[j+1] - ofill[j]); 1193c0ab637bSBarry Smith maxcnt = PetscMax(maxcnt,ocols[cnt]+cols[cnt]); 1194ce308e1dSBarry Smith cnt++; 1195ce308e1dSBarry Smith } 1196ce308e1dSBarry Smith } 1197ce308e1dSBarry Smith /* coupling with process to the right */ 1198ce308e1dSBarry Smith for (i=nx-s; i<nx; i++) { 1199ce308e1dSBarry Smith for (j=0; j<nc; j++) { 1200ce308e1dSBarry Smith ocols[cnt] = ((rank == (size-1)) ? 0 : (i - nx + s + 1)*(ofill[j+1] - ofill[j])); 12010acb5bebSBarry Smith cols[cnt] = dfill[j+1] - dfill[j] + (s + nx - i - 1)*(ofill[j+1] - ofill[j]); 1202c0ab637bSBarry Smith maxcnt = PetscMax(maxcnt,ocols[cnt]+cols[cnt]); 1203ce308e1dSBarry Smith cnt++; 1204ce308e1dSBarry Smith } 1205ce308e1dSBarry Smith } 1206ce308e1dSBarry Smith 1207ce308e1dSBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,cols);CHKERRQ(ierr); 1208ce308e1dSBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,cols,0,ocols);CHKERRQ(ierr); 1209ce308e1dSBarry Smith ierr = PetscFree2(cols,ocols);CHKERRQ(ierr); 1210ce308e1dSBarry Smith 1211ce308e1dSBarry Smith ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 1212ce308e1dSBarry Smith ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 1213ce308e1dSBarry Smith 1214ce308e1dSBarry Smith /* 1215ce308e1dSBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 1216ce308e1dSBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 1217ce308e1dSBarry Smith PETSc ordering. 1218ce308e1dSBarry Smith */ 1219ce308e1dSBarry Smith if (!da->prealloc_only) { 1220c0ab637bSBarry Smith ierr = PetscCalloc2(maxcnt,&values,maxcnt,&cols);CHKERRQ(ierr); 1221ce308e1dSBarry Smith 1222ce308e1dSBarry Smith row = xs*nc; 1223ce308e1dSBarry Smith /* coupling with process to the left */ 1224ce308e1dSBarry Smith for (i=xs; i<xs+s; i++) { 1225ce308e1dSBarry Smith for (j=0; j<nc; j++) { 1226ce308e1dSBarry Smith cnt = 0; 1227ce308e1dSBarry Smith if (rank) { 1228ce308e1dSBarry Smith for (l=0; l<s; l++) { 1229ce308e1dSBarry Smith for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k]; 1230ce308e1dSBarry Smith } 1231ce308e1dSBarry Smith } 12320acb5bebSBarry Smith if (dfill) { 12330acb5bebSBarry Smith for (k=dfill[j]; k<dfill[j+1]; k++) { 12340acb5bebSBarry Smith cols[cnt++] = i*nc + dfill[k]; 12350acb5bebSBarry Smith } 12360acb5bebSBarry Smith } else { 1237ce308e1dSBarry Smith for (k=0; k<nc; k++) { 1238ce308e1dSBarry Smith cols[cnt++] = i*nc + k; 1239ce308e1dSBarry Smith } 12400acb5bebSBarry Smith } 1241ce308e1dSBarry Smith for (l=0; l<s; l++) { 1242ce308e1dSBarry Smith for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k]; 1243ce308e1dSBarry Smith } 1244ce308e1dSBarry Smith ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 1245ce308e1dSBarry Smith row++; 1246ce308e1dSBarry Smith } 1247ce308e1dSBarry Smith } 1248ce308e1dSBarry Smith for (i=xs+s; i<xs+nx-s; i++) { 1249ce308e1dSBarry Smith for (j=0; j<nc; j++) { 1250ce308e1dSBarry Smith cnt = 0; 1251ce308e1dSBarry Smith for (l=0; l<s; l++) { 1252ce308e1dSBarry Smith for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k]; 1253ce308e1dSBarry Smith } 12540acb5bebSBarry Smith if (dfill) { 12550acb5bebSBarry Smith for (k=dfill[j]; k<dfill[j+1]; k++) { 12560acb5bebSBarry Smith cols[cnt++] = i*nc + dfill[k]; 12570acb5bebSBarry Smith } 12580acb5bebSBarry Smith } else { 1259ce308e1dSBarry Smith for (k=0; k<nc; k++) { 1260ce308e1dSBarry Smith cols[cnt++] = i*nc + k; 1261ce308e1dSBarry Smith } 12620acb5bebSBarry Smith } 1263ce308e1dSBarry Smith for (l=0; l<s; l++) { 1264ce308e1dSBarry Smith for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k]; 1265ce308e1dSBarry Smith } 1266ce308e1dSBarry Smith ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 1267ce308e1dSBarry Smith row++; 1268ce308e1dSBarry Smith } 1269ce308e1dSBarry Smith } 1270ce308e1dSBarry Smith /* coupling with process to the right */ 1271ce308e1dSBarry Smith for (i=xs+nx-s; i<xs+nx; i++) { 1272ce308e1dSBarry Smith for (j=0; j<nc; j++) { 1273ce308e1dSBarry Smith cnt = 0; 1274ce308e1dSBarry Smith for (l=0; l<s; l++) { 1275ce308e1dSBarry Smith for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k]; 1276ce308e1dSBarry Smith } 12770acb5bebSBarry Smith if (dfill) { 12780acb5bebSBarry Smith for (k=dfill[j]; k<dfill[j+1]; k++) { 12790acb5bebSBarry Smith cols[cnt++] = i*nc + dfill[k]; 12800acb5bebSBarry Smith } 12810acb5bebSBarry Smith } else { 1282ce308e1dSBarry Smith for (k=0; k<nc; k++) { 1283ce308e1dSBarry Smith cols[cnt++] = i*nc + k; 1284ce308e1dSBarry Smith } 12850acb5bebSBarry Smith } 1286ce308e1dSBarry Smith if (rank < size-1) { 1287ce308e1dSBarry Smith for (l=0; l<s; l++) { 1288ce308e1dSBarry Smith for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k]; 1289ce308e1dSBarry Smith } 1290ce308e1dSBarry Smith } 1291ce308e1dSBarry Smith ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 1292ce308e1dSBarry Smith row++; 1293ce308e1dSBarry Smith } 1294ce308e1dSBarry Smith } 1295c0ab637bSBarry Smith ierr = PetscFree2(values,cols);CHKERRQ(ierr); 1296ce308e1dSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1297ce308e1dSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1298189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1299ce308e1dSBarry Smith } 1300ce308e1dSBarry Smith PetscFunctionReturn(0); 1301ce308e1dSBarry Smith } 1302ce308e1dSBarry Smith 1303ce308e1dSBarry Smith /* ---------------------------------------------------------------------------------*/ 1304ce308e1dSBarry Smith 1305950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM da,Mat J) 130647c6ae99SBarry Smith { 130747c6ae99SBarry Smith PetscErrorCode ierr; 130847c6ae99SBarry Smith PetscInt xs,nx,i,i1,slot,gxs,gnx; 13090298fd71SBarry Smith PetscInt m,dim,s,*cols = NULL,nc,*rows = NULL,col,cnt,l; 131047c6ae99SBarry Smith PetscInt istart,iend; 131147c6ae99SBarry Smith PetscScalar *values; 1312bff4a2f0SMatthew G. Knepley DMBoundaryType bx; 131345b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 131447c6ae99SBarry Smith 131547c6ae99SBarry Smith PetscFunctionBegin; 131647c6ae99SBarry Smith /* 131747c6ae99SBarry Smith nc - number of components per grid point 131847c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 131947c6ae99SBarry Smith 132047c6ae99SBarry Smith */ 13211321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr); 132247c6ae99SBarry Smith col = 2*s + 1; 132347c6ae99SBarry Smith 1324aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr); 1325aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr); 132647c6ae99SBarry Smith 1327f73d5cc4SBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 132847c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr); 132947c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr); 133047c6ae99SBarry Smith 13311411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 1332784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 133347c6ae99SBarry Smith 133447c6ae99SBarry Smith /* 133547c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 133647c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 133747c6ae99SBarry Smith PETSc ordering. 133847c6ae99SBarry Smith */ 1339fcfd50ebSBarry Smith if (!da->prealloc_only) { 1340dcca6d9dSJed Brown ierr = PetscMalloc2(nc,&rows,col*nc*nc,&cols);CHKERRQ(ierr); 13411795a4d1SJed Brown ierr = PetscCalloc1(col*nc*nc,&values);CHKERRQ(ierr); 134247c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 134347c6ae99SBarry Smith istart = PetscMax(-s,gxs - i); 134447c6ae99SBarry Smith iend = PetscMin(s,gxs + gnx - i - 1); 134547c6ae99SBarry Smith slot = i - gxs; 134647c6ae99SBarry Smith 134747c6ae99SBarry Smith cnt = 0; 134847c6ae99SBarry Smith for (l=0; l<nc; l++) { 134947c6ae99SBarry Smith for (i1=istart; i1<iend+1; i1++) { 135047c6ae99SBarry Smith cols[cnt++] = l + nc*(slot + i1); 135147c6ae99SBarry Smith } 135247c6ae99SBarry Smith rows[l] = l + nc*(slot); 135347c6ae99SBarry Smith } 135447c6ae99SBarry Smith ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 135547c6ae99SBarry Smith } 135647c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 135747c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 135847c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1359189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 136047c6ae99SBarry Smith ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 1361ce308e1dSBarry Smith } 136247c6ae99SBarry Smith PetscFunctionReturn(0); 136347c6ae99SBarry Smith } 136447c6ae99SBarry Smith 1365950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM da,Mat J) 136647c6ae99SBarry Smith { 136747c6ae99SBarry Smith PetscErrorCode ierr; 136847c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 136947c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz; 137047c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,ii,jj; 137147c6ae99SBarry Smith MPI_Comm comm; 137247c6ae99SBarry Smith PetscScalar *values; 1373bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by; 1374aa219208SBarry Smith DMDAStencilType st; 137545b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 137647c6ae99SBarry Smith 137747c6ae99SBarry Smith PetscFunctionBegin; 137847c6ae99SBarry Smith /* 137947c6ae99SBarry Smith nc - number of components per grid point 138047c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 138147c6ae99SBarry Smith */ 13821321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 138347c6ae99SBarry Smith col = 2*s + 1; 138447c6ae99SBarry Smith 1385aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 1386aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 138747c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 138847c6ae99SBarry Smith 1389785e854fSJed Brown ierr = PetscMalloc1(col*col*nc*nc,&cols);CHKERRQ(ierr); 139047c6ae99SBarry Smith 13911411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 139247c6ae99SBarry Smith 139347c6ae99SBarry Smith /* determine the matrix preallocation information */ 139447c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr); 139547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1396bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1397bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 139847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1399bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1400bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 140147c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 140247c6ae99SBarry Smith 140347c6ae99SBarry Smith /* Find block columns in block row */ 140447c6ae99SBarry Smith cnt = 0; 140547c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 140647c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1407aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */ 140847c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 140947c6ae99SBarry Smith } 141047c6ae99SBarry Smith } 141147c6ae99SBarry Smith } 1412d6e23781SBarry Smith ierr = MatPreallocateSetLocalBlock(ltog,1,&slot,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 141347c6ae99SBarry Smith } 141447c6ae99SBarry Smith } 141547c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 141647c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 141747c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 141847c6ae99SBarry Smith 1419784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 142047c6ae99SBarry Smith 142147c6ae99SBarry Smith /* 142247c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 142347c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 142447c6ae99SBarry Smith PETSc ordering. 142547c6ae99SBarry Smith */ 1426fcfd50ebSBarry Smith if (!da->prealloc_only) { 14271795a4d1SJed Brown ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr); 142847c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1429bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1430bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 143147c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1432bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1433bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 143447c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 143547c6ae99SBarry Smith cnt = 0; 143647c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 143747c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1438aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */ 143947c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 144047c6ae99SBarry Smith } 144147c6ae99SBarry Smith } 144247c6ae99SBarry Smith } 144347c6ae99SBarry Smith ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 144447c6ae99SBarry Smith } 144547c6ae99SBarry Smith } 144647c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 144747c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 144847c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1449189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 145047c6ae99SBarry Smith } 145147c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 145247c6ae99SBarry Smith PetscFunctionReturn(0); 145347c6ae99SBarry Smith } 145447c6ae99SBarry Smith 1455950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM da,Mat J) 145647c6ae99SBarry Smith { 145747c6ae99SBarry Smith PetscErrorCode ierr; 145847c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 145947c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz; 146047c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 146147c6ae99SBarry Smith MPI_Comm comm; 146247c6ae99SBarry Smith PetscScalar *values; 1463bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by,bz; 1464aa219208SBarry Smith DMDAStencilType st; 146545b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 146647c6ae99SBarry Smith 146747c6ae99SBarry Smith PetscFunctionBegin; 146847c6ae99SBarry Smith /* 146947c6ae99SBarry Smith nc - number of components per grid point 147047c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 147147c6ae99SBarry Smith 147247c6ae99SBarry Smith */ 14731321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 147447c6ae99SBarry Smith col = 2*s + 1; 147547c6ae99SBarry Smith 1476aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1477aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 147847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 147947c6ae99SBarry Smith 1480785e854fSJed Brown ierr = PetscMalloc1(col*col*col,&cols);CHKERRQ(ierr); 148147c6ae99SBarry Smith 14821411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 148347c6ae99SBarry Smith 148447c6ae99SBarry Smith /* determine the matrix preallocation information */ 148547c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr); 148647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1487bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1488bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 148947c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1490bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1491bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 149247c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1493bff4a2f0SMatthew G. Knepley kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 1494bff4a2f0SMatthew G. Knepley kend = (bz == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 149547c6ae99SBarry Smith 149647c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 149747c6ae99SBarry Smith 149847c6ae99SBarry Smith /* Find block columns in block row */ 149947c6ae99SBarry Smith cnt = 0; 150047c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 150147c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 150247c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1503aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 150447c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 150547c6ae99SBarry Smith } 150647c6ae99SBarry Smith } 150747c6ae99SBarry Smith } 150847c6ae99SBarry Smith } 1509d6e23781SBarry Smith ierr = MatPreallocateSetLocalBlock(ltog,1,&slot,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 151047c6ae99SBarry Smith } 151147c6ae99SBarry Smith } 151247c6ae99SBarry Smith } 151347c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 151447c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 151547c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 151647c6ae99SBarry Smith 1517784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 151847c6ae99SBarry Smith 151947c6ae99SBarry Smith /* 152047c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 152147c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 152247c6ae99SBarry Smith PETSc ordering. 152347c6ae99SBarry Smith */ 1524fcfd50ebSBarry Smith if (!da->prealloc_only) { 15251795a4d1SJed Brown ierr = PetscCalloc1(col*col*col*nc*nc,&values);CHKERRQ(ierr); 152647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1527bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1528bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 152947c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1530bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1531bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 153247c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1533bff4a2f0SMatthew G. Knepley kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 1534bff4a2f0SMatthew G. Knepley kend = (bz == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 153547c6ae99SBarry Smith 153647c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 153747c6ae99SBarry Smith 153847c6ae99SBarry Smith cnt = 0; 153947c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 154047c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 154147c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1542aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 154347c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 154447c6ae99SBarry Smith } 154547c6ae99SBarry Smith } 154647c6ae99SBarry Smith } 154747c6ae99SBarry Smith } 154847c6ae99SBarry Smith ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 154947c6ae99SBarry Smith } 155047c6ae99SBarry Smith } 155147c6ae99SBarry Smith } 155247c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 155347c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 155447c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1555189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 155647c6ae99SBarry Smith } 155747c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 155847c6ae99SBarry Smith PetscFunctionReturn(0); 155947c6ae99SBarry Smith } 156047c6ae99SBarry Smith 156147c6ae99SBarry Smith /* 156247c6ae99SBarry Smith This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to 156347c6ae99SBarry Smith identify in the local ordering with periodic domain. 156447c6ae99SBarry Smith */ 156547c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[]) 156647c6ae99SBarry Smith { 156747c6ae99SBarry Smith PetscErrorCode ierr; 156847c6ae99SBarry Smith PetscInt i,n; 156947c6ae99SBarry Smith 157047c6ae99SBarry Smith PetscFunctionBegin; 1571d6e23781SBarry Smith ierr = ISLocalToGlobalMappingApplyBlock(ltog,1,row,row);CHKERRQ(ierr); 1572d6e23781SBarry Smith ierr = ISLocalToGlobalMappingApplyBlock(ltog,*cnt,col,col);CHKERRQ(ierr); 157347c6ae99SBarry Smith for (i=0,n=0; i<*cnt; i++) { 157447c6ae99SBarry Smith if (col[i] >= *row) col[n++] = col[i]; 157547c6ae99SBarry Smith } 157647c6ae99SBarry Smith *cnt = n; 157747c6ae99SBarry Smith PetscFunctionReturn(0); 157847c6ae99SBarry Smith } 157947c6ae99SBarry Smith 1580950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM da,Mat J) 158147c6ae99SBarry Smith { 158247c6ae99SBarry Smith PetscErrorCode ierr; 158347c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 158447c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz; 158547c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,ii,jj; 158647c6ae99SBarry Smith MPI_Comm comm; 158747c6ae99SBarry Smith PetscScalar *values; 1588bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by; 1589aa219208SBarry Smith DMDAStencilType st; 159045b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 159147c6ae99SBarry Smith 159247c6ae99SBarry Smith PetscFunctionBegin; 159347c6ae99SBarry Smith /* 159447c6ae99SBarry Smith nc - number of components per grid point 159547c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 159647c6ae99SBarry Smith */ 15971321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr); 159847c6ae99SBarry Smith col = 2*s + 1; 159947c6ae99SBarry Smith 1600aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr); 1601aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr); 160247c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 160347c6ae99SBarry Smith 1604785e854fSJed Brown ierr = PetscMalloc1(col*col*nc*nc,&cols);CHKERRQ(ierr); 160547c6ae99SBarry Smith 16061411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 160747c6ae99SBarry Smith 160847c6ae99SBarry Smith /* determine the matrix preallocation information */ 1609eabe889fSLisandro Dalcin ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr); 161047c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1611bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1612bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 161347c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1614bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1615bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 161647c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 161747c6ae99SBarry Smith 161847c6ae99SBarry Smith /* Find block columns in block row */ 161947c6ae99SBarry Smith cnt = 0; 162047c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 162147c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1622aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { 162347c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 162447c6ae99SBarry Smith } 162547c6ae99SBarry Smith } 162647c6ae99SBarry Smith } 162745b6f7e9SBarry Smith ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr); 1628d6e23781SBarry Smith ierr = MatPreallocateSymmetricSetBlock(slot,cnt,cols,dnz,onz);CHKERRQ(ierr); 162947c6ae99SBarry Smith } 163047c6ae99SBarry Smith } 163147c6ae99SBarry Smith ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 163247c6ae99SBarry Smith ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 163347c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 163447c6ae99SBarry Smith 1635784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 163647c6ae99SBarry Smith 163747c6ae99SBarry Smith /* 163847c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 163947c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 164047c6ae99SBarry Smith PETSc ordering. 164147c6ae99SBarry Smith */ 1642fcfd50ebSBarry Smith if (!da->prealloc_only) { 16431795a4d1SJed Brown ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr); 164447c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1645bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1646bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 164747c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1648bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1649bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 165047c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys); 165147c6ae99SBarry Smith 165247c6ae99SBarry Smith /* Find block columns in block row */ 165347c6ae99SBarry Smith cnt = 0; 165447c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 165547c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 1656aa219208SBarry Smith if (st == DMDA_STENCIL_BOX || !ii || !jj) { 165747c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj; 165847c6ae99SBarry Smith } 165947c6ae99SBarry Smith } 166047c6ae99SBarry Smith } 166145b6f7e9SBarry Smith ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr); 166247c6ae99SBarry Smith ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 166347c6ae99SBarry Smith } 166447c6ae99SBarry Smith } 166547c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 166647c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166747c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1668189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 166947c6ae99SBarry Smith } 167047c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 167147c6ae99SBarry Smith PetscFunctionReturn(0); 167247c6ae99SBarry Smith } 167347c6ae99SBarry Smith 1674950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM da,Mat J) 167547c6ae99SBarry Smith { 167647c6ae99SBarry Smith PetscErrorCode ierr; 167747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 167847c6ae99SBarry Smith PetscInt m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz; 167947c6ae99SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk; 168047c6ae99SBarry Smith MPI_Comm comm; 168147c6ae99SBarry Smith PetscScalar *values; 1682bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by,bz; 1683aa219208SBarry Smith DMDAStencilType st; 168445b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 168547c6ae99SBarry Smith 168647c6ae99SBarry Smith PetscFunctionBegin; 168747c6ae99SBarry Smith /* 168847c6ae99SBarry Smith nc - number of components per grid point 168947c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 169047c6ae99SBarry Smith */ 16911321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 169247c6ae99SBarry Smith col = 2*s + 1; 169347c6ae99SBarry Smith 1694aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1695aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 169647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 169747c6ae99SBarry Smith 169847c6ae99SBarry Smith /* create the matrix */ 1699785e854fSJed Brown ierr = PetscMalloc1(col*col*col,&cols);CHKERRQ(ierr); 170047c6ae99SBarry Smith 17011411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 170247c6ae99SBarry Smith 170347c6ae99SBarry Smith /* determine the matrix preallocation information */ 1704eabe889fSLisandro Dalcin ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr); 170547c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1706bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1707bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 170847c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1709bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1710bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 171147c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1712bff4a2f0SMatthew G. Knepley kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 1713bff4a2f0SMatthew G. Knepley kend = (bz == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 171447c6ae99SBarry Smith 171547c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 171647c6ae99SBarry Smith 171747c6ae99SBarry Smith /* Find block columns in block row */ 171847c6ae99SBarry Smith cnt = 0; 171947c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 172047c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 172147c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1722aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) { 172347c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 172447c6ae99SBarry Smith } 172547c6ae99SBarry Smith } 172647c6ae99SBarry Smith } 172747c6ae99SBarry Smith } 172845b6f7e9SBarry Smith ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr); 1729d6e23781SBarry Smith ierr = MatPreallocateSymmetricSetBlock(slot,cnt,cols,dnz,onz);CHKERRQ(ierr); 173047c6ae99SBarry Smith } 173147c6ae99SBarry Smith } 173247c6ae99SBarry Smith } 173347c6ae99SBarry Smith ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr); 173447c6ae99SBarry Smith ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr); 173547c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 173647c6ae99SBarry Smith 1737784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 173847c6ae99SBarry Smith 173947c6ae99SBarry Smith /* 174047c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 174147c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 174247c6ae99SBarry Smith PETSc ordering. 174347c6ae99SBarry Smith */ 1744fcfd50ebSBarry Smith if (!da->prealloc_only) { 17451795a4d1SJed Brown ierr = PetscCalloc1(col*col*col*nc*nc,&values);CHKERRQ(ierr); 174647c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1747bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1748bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 174947c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1750bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1751bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 175247c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1753bff4a2f0SMatthew G. Knepley kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 1754bff4a2f0SMatthew G. Knepley kend = (bz == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 175547c6ae99SBarry Smith 175647c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 175747c6ae99SBarry Smith 175847c6ae99SBarry Smith cnt = 0; 175947c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 176047c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 176147c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 1762aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) { 176347c6ae99SBarry Smith cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; 176447c6ae99SBarry Smith } 176547c6ae99SBarry Smith } 176647c6ae99SBarry Smith } 176747c6ae99SBarry Smith } 176845b6f7e9SBarry Smith ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr); 176947c6ae99SBarry Smith ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 177047c6ae99SBarry Smith } 177147c6ae99SBarry Smith } 177247c6ae99SBarry Smith } 177347c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 177447c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 177547c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1776189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 177747c6ae99SBarry Smith } 177847c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 177947c6ae99SBarry Smith PetscFunctionReturn(0); 178047c6ae99SBarry Smith } 178147c6ae99SBarry Smith 178247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/ 178347c6ae99SBarry Smith 1784950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J) 178547c6ae99SBarry Smith { 178647c6ae99SBarry Smith PetscErrorCode ierr; 178747c6ae99SBarry Smith PetscInt xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny; 1788c0ab637bSBarry Smith PetscInt m,n,dim,s,*cols,k,nc,row,col,cnt, maxcnt = 0,l,p,*dnz,*onz; 1789c1154cd5SBarry Smith PetscInt istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk,M,N,P; 179047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 179147c6ae99SBarry Smith PetscInt ifill_col,*dfill = dd->dfill,*ofill = dd->ofill; 179247c6ae99SBarry Smith MPI_Comm comm; 179347c6ae99SBarry Smith PetscScalar *values; 1794bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by,bz; 179545b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 1796aa219208SBarry Smith DMDAStencilType st; 1797c1154cd5SBarry Smith PetscBool removedups = PETSC_FALSE; 179847c6ae99SBarry Smith 179947c6ae99SBarry Smith PetscFunctionBegin; 180047c6ae99SBarry Smith /* 180147c6ae99SBarry Smith nc - number of components per grid point 180247c6ae99SBarry Smith col - number of colors needed in one direction for single component problem 180347c6ae99SBarry Smith 180447c6ae99SBarry Smith */ 1805c1154cd5SBarry Smith ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr); 180647c6ae99SBarry Smith col = 2*s + 1; 1807bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\ 180847c6ae99SBarry Smith by 2*stencil_width + 1\n"); 1809bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\ 181047c6ae99SBarry Smith by 2*stencil_width + 1\n"); 1811bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\ 181247c6ae99SBarry Smith by 2*stencil_width + 1\n"); 181347c6ae99SBarry Smith 1814c1154cd5SBarry Smith /* 1815c1154cd5SBarry Smith With one processor in periodic domains in a skinny dimension the code will label nonzero columns multiple times 1816c1154cd5SBarry Smith because of "wrapping" around the end of the domain hitting an entry already counted in the other direction. 1817c1154cd5SBarry Smith */ 1818c1154cd5SBarry Smith if (M == 1 && 2*s >= m) removedups = PETSC_TRUE; 1819c1154cd5SBarry Smith if (N == 1 && 2*s >= n) removedups = PETSC_TRUE; 1820c1154cd5SBarry Smith if (P == 1 && 2*s >= p) removedups = PETSC_TRUE; 1821c1154cd5SBarry Smith 1822aa219208SBarry Smith ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr); 1823aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr); 182447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 182547c6ae99SBarry Smith 1826785e854fSJed Brown ierr = PetscMalloc1(col*col*col*nc,&cols);CHKERRQ(ierr); 18271411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<og);CHKERRQ(ierr); 182847c6ae99SBarry Smith 182947c6ae99SBarry Smith /* determine the matrix preallocation information */ 183047c6ae99SBarry Smith ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr); 183147c6ae99SBarry Smith 183206ca8cadSBarry Smith ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr); 183347c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1834bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1835bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 183647c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1837bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1838bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 183947c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1840bff4a2f0SMatthew G. Knepley kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 1841bff4a2f0SMatthew G. Knepley kend = (bz == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 184247c6ae99SBarry Smith 184347c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 184447c6ae99SBarry Smith 184547c6ae99SBarry Smith for (l=0; l<nc; l++) { 184647c6ae99SBarry Smith cnt = 0; 184747c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 184847c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 184947c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 185047c6ae99SBarry Smith if (ii || jj || kk) { 1851aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 18528865f1eaSKarl Rupp for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 185347c6ae99SBarry Smith } 185447c6ae99SBarry Smith } else { 185547c6ae99SBarry Smith if (dfill) { 18568865f1eaSKarl Rupp for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 185747c6ae99SBarry Smith } else { 18588865f1eaSKarl Rupp for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk); 185947c6ae99SBarry Smith } 186047c6ae99SBarry Smith } 186147c6ae99SBarry Smith } 186247c6ae99SBarry Smith } 186347c6ae99SBarry Smith } 186447c6ae99SBarry Smith row = l + nc*(slot); 1865c0ab637bSBarry Smith maxcnt = PetscMax(maxcnt,cnt); 1866c1154cd5SBarry Smith if (removedups) { 1867c1154cd5SBarry Smith ierr = MatPreallocateSetLocalRemoveDups(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 1868c1154cd5SBarry Smith } else { 1869784ac674SJed Brown ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr); 187047c6ae99SBarry Smith } 187147c6ae99SBarry Smith } 187247c6ae99SBarry Smith } 187347c6ae99SBarry Smith } 1874c1154cd5SBarry Smith } 187547c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 187647c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 187747c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1878784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr); 187947c6ae99SBarry Smith 188047c6ae99SBarry Smith /* 188147c6ae99SBarry Smith For each node in the grid: we get the neighbors in the local (on processor ordering 188247c6ae99SBarry Smith that includes the ghost points) then MatSetValuesLocal() maps those indices to the global 188347c6ae99SBarry Smith PETSc ordering. 188447c6ae99SBarry Smith */ 1885fcfd50ebSBarry Smith if (!da->prealloc_only) { 1886c0ab637bSBarry Smith ierr = PetscCalloc1(maxcnt,&values);CHKERRQ(ierr); 188747c6ae99SBarry Smith for (i=xs; i<xs+nx; i++) { 1888bff4a2f0SMatthew G. Knepley istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i)); 1889bff4a2f0SMatthew G. Knepley iend = (bx == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,m-i-1)); 189047c6ae99SBarry Smith for (j=ys; j<ys+ny; j++) { 1891bff4a2f0SMatthew G. Knepley jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j)); 1892bff4a2f0SMatthew G. Knepley jend = (by == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,n-j-1)); 189347c6ae99SBarry Smith for (k=zs; k<zs+nz; k++) { 1894bff4a2f0SMatthew G. Knepley kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k)); 1895bff4a2f0SMatthew G. Knepley kend = (bz == DM_BOUNDARY_PERIODIC) ? s : (PetscMin(s,p-k-1)); 189647c6ae99SBarry Smith 189747c6ae99SBarry Smith slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs); 189847c6ae99SBarry Smith 189947c6ae99SBarry Smith for (l=0; l<nc; l++) { 190047c6ae99SBarry Smith cnt = 0; 190147c6ae99SBarry Smith for (ii=istart; ii<iend+1; ii++) { 190247c6ae99SBarry Smith for (jj=jstart; jj<jend+1; jj++) { 190347c6ae99SBarry Smith for (kk=kstart; kk<kend+1; kk++) { 190447c6ae99SBarry Smith if (ii || jj || kk) { 1905aa219208SBarry Smith if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/ 19068865f1eaSKarl Rupp for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 190747c6ae99SBarry Smith } 190847c6ae99SBarry Smith } else { 190947c6ae99SBarry Smith if (dfill) { 19108865f1eaSKarl Rupp for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk); 191147c6ae99SBarry Smith } else { 19128865f1eaSKarl Rupp for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk); 191347c6ae99SBarry Smith } 191447c6ae99SBarry Smith } 191547c6ae99SBarry Smith } 191647c6ae99SBarry Smith } 191747c6ae99SBarry Smith } 191847c6ae99SBarry Smith row = l + nc*(slot); 191947c6ae99SBarry Smith ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr); 192047c6ae99SBarry Smith } 192147c6ae99SBarry Smith } 192247c6ae99SBarry Smith } 192347c6ae99SBarry Smith } 192447c6ae99SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 192547c6ae99SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 192647c6ae99SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1927189e4007SBarry Smith ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 192847c6ae99SBarry Smith } 192947c6ae99SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 193047c6ae99SBarry Smith PetscFunctionReturn(0); 193147c6ae99SBarry Smith } 1932