17d0a6c19SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular 1d arrays in parallel. 447c6ae99SBarry Smith This file was created by Peter Mell 6/30/95 547c6ae99SBarry Smith */ 647c6ae99SBarry Smith 7af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 847c6ae99SBarry Smith 99804daf3SBarry Smith #include <petscdraw.h> 10e0877f53SBarry Smith static PetscErrorCode DMView_DA_1d(DM da,PetscViewer viewer) 1147c6ae99SBarry Smith { 1247c6ae99SBarry Smith PetscErrorCode ierr; 1347c6ae99SBarry Smith PetscMPIInt rank; 14c9493c35SStefano Zampini PetscBool iascii,isdraw,isglvis,isbinary; 1547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 169a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 179a42bb27SBarry Smith PetscBool ismatlab; 189a42bb27SBarry Smith #endif 1947c6ae99SBarry Smith 2047c6ae99SBarry Smith PetscFunctionBegin; 21ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr); 2247c6ae99SBarry Smith 23251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 24251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 25c9493c35SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr); 26251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 279a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 28251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr); 299a42bb27SBarry Smith #endif 3047c6ae99SBarry Smith if (iascii) { 3147c6ae99SBarry Smith PetscViewerFormat format; 3247c6ae99SBarry Smith 3347c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 3476a8abe0SBarry Smith if (format == PETSC_VIEWER_LOAD_BALANCE) { 3576a8abe0SBarry Smith PetscInt i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal; 3676a8abe0SBarry Smith DMDALocalInfo info; 3776a8abe0SBarry Smith PetscMPIInt size; 3876a8abe0SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr); 3976a8abe0SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 4076a8abe0SBarry Smith nzlocal = info.xm; 4176a8abe0SBarry Smith ierr = PetscMalloc1(size,&nz);CHKERRQ(ierr); 4276a8abe0SBarry Smith ierr = MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 4376a8abe0SBarry Smith for (i=0; i<(PetscInt)size; i++) { 4476a8abe0SBarry Smith nmax = PetscMax(nmax,nz[i]); 4576a8abe0SBarry Smith nmin = PetscMin(nmin,nz[i]); 4676a8abe0SBarry Smith navg += nz[i]; 4776a8abe0SBarry Smith } 4876a8abe0SBarry Smith ierr = PetscFree(nz);CHKERRQ(ierr); 4976a8abe0SBarry Smith navg = navg/size; 5076a8abe0SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Load Balance - Grid Points: Min %D avg %D max %D\n",nmin,navg,nmax);CHKERRQ(ierr); 5176a8abe0SBarry Smith PetscFunctionReturn(0); 5276a8abe0SBarry Smith } 538135c375SStefano Zampini if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL && format != PETSC_VIEWER_ASCII_GLVIS) { 54aa219208SBarry Smith DMDALocalInfo info; 55aa219208SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 561575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 5747c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D m %D w %D s %D\n",rank,dd->M,dd->m,dd->w,dd->s);CHKERRQ(ierr); 5847c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D\n",info.xs,info.xs+info.xm);CHKERRQ(ierr); 5947c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 601575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 618135c375SStefano Zampini } else if (format == PETSC_VIEWER_ASCII_GLVIS) { 628135c375SStefano Zampini ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr); 639a42bb27SBarry Smith } else { 649a42bb27SBarry Smith ierr = DMView_DA_VTK(da, viewer);CHKERRQ(ierr); 6547c6ae99SBarry Smith } 6647c6ae99SBarry Smith } else if (isdraw) { 6747c6ae99SBarry Smith PetscDraw draw; 6847c6ae99SBarry Smith double ymin = -1,ymax = 1,xmin = -1,xmax = dd->M,x; 6947c6ae99SBarry Smith PetscInt base; 7047c6ae99SBarry Smith char node[10]; 7147c6ae99SBarry Smith PetscBool isnull; 7247c6ae99SBarry Smith 7347c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 745b399a63SLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 755b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 765b399a63SLisandro Dalcin 775b399a63SLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 785b399a63SLisandro Dalcin ierr = PetscDrawClear(draw);CHKERRQ(ierr); 7947c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 8047c6ae99SBarry Smith 815b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 8247c6ae99SBarry Smith /* first processor draws all node lines */ 8347c6ae99SBarry Smith if (!rank) { 8447c6ae99SBarry Smith PetscInt xmin_tmp; 8547c6ae99SBarry Smith ymin = 0.0; ymax = 0.3; 8647c6ae99SBarry Smith for (xmin_tmp=0; xmin_tmp < dd->M; xmin_tmp++) { 8747c6ae99SBarry Smith ierr = PetscDrawLine(draw,(double)xmin_tmp,ymin,(double)xmin_tmp,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 8847c6ae99SBarry Smith } 8947c6ae99SBarry Smith xmin = 0.0; xmax = dd->M - 1; 9047c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 9147c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 9247c6ae99SBarry Smith } 935b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 945b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 9547c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 9647c6ae99SBarry Smith 975b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 9847c6ae99SBarry Smith /* draw my box */ 9947c6ae99SBarry Smith ymin = 0; ymax = 0.3; xmin = dd->xs / dd->w; xmax = (dd->xe / dd->w) - 1; 10047c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 10147c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10247c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10347c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10447c6ae99SBarry Smith /* Put in index numbers */ 10547c6ae99SBarry Smith base = dd->base / dd->w; 10647c6ae99SBarry Smith for (x=xmin; x<=xmax; x++) { 107832b7cebSLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr); 10847c6ae99SBarry Smith ierr = PetscDrawString(draw,x,ymin,PETSC_DRAW_RED,node);CHKERRQ(ierr); 10947c6ae99SBarry Smith } 1105b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1115b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 11247c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 113832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 114c9493c35SStefano Zampini } else if (isglvis) { 115c9493c35SStefano Zampini ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr); 1169a42bb27SBarry Smith } else if (isbinary) { 1179a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1189a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1199a42bb27SBarry Smith } else if (ismatlab) { 1209a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1219a42bb27SBarry Smith #endif 12211aeaf0aSBarry Smith } 12347c6ae99SBarry Smith PetscFunctionReturn(0); 12447c6ae99SBarry Smith } 12547c6ae99SBarry Smith 12647c6ae99SBarry Smith 1277087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_1D(DM da) 12847c6ae99SBarry Smith { 12947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 13047c6ae99SBarry Smith const PetscInt M = dd->M; 13147c6ae99SBarry Smith const PetscInt dof = dd->w; 13247c6ae99SBarry Smith const PetscInt s = dd->s; 13345b6f7e9SBarry Smith const PetscInt sDist = s; /* stencil distance in points */ 13447c6ae99SBarry Smith const PetscInt *lx = dd->lx; 135bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 13647c6ae99SBarry Smith MPI_Comm comm; 13747c6ae99SBarry Smith Vec local, global; 138bd1fc5aeSBarry Smith VecScatter gtol; 13947c6ae99SBarry Smith IS to, from; 14047c6ae99SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE; 14147c6ae99SBarry Smith PetscMPIInt rank, size; 142bd1fc5aeSBarry Smith PetscInt i,*idx,nn,left,xs,xe,x,Xs,Xe,start,m,IXs,IXe; 14347c6ae99SBarry Smith PetscErrorCode ierr; 14447c6ae99SBarry Smith 14547c6ae99SBarry Smith PetscFunctionBegin; 14647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 14747c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 14847c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 14947c6ae99SBarry Smith 1507d310018SBarry Smith dd->p = 1; 1517d310018SBarry Smith dd->n = 1; 15247c6ae99SBarry Smith dd->m = size; 15347c6ae99SBarry Smith m = dd->m; 15447c6ae99SBarry Smith 15547c6ae99SBarry Smith if (s > 0) { 15647c6ae99SBarry Smith /* if not communicating data then should be ok to have nothing on some processes */ 15747c6ae99SBarry Smith if (M < m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"More processes than data points! %D %D",m,M); 158a66092f1SBarry Smith if ((M-1) < s && size > 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array is too small for stencil! %D %D",M-1,s); 15947c6ae99SBarry Smith } 16047c6ae99SBarry Smith 16147c6ae99SBarry Smith /* 16247c6ae99SBarry Smith Determine locally owned region 16347c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 16447c6ae99SBarry Smith */ 16547c6ae99SBarry Smith if (!lx) { 166785e854fSJed Brown ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr); 167c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_partition_blockcomm",&flg1,NULL);CHKERRQ(ierr); 168c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_partition_nodes_at_end",&flg2,NULL);CHKERRQ(ierr); 16947c6ae99SBarry Smith if (flg1) { /* Block Comm type Distribution */ 17047c6ae99SBarry Smith xs = rank*M/m; 17147c6ae99SBarry Smith x = (rank + 1)*M/m - xs; 17247c6ae99SBarry Smith } else if (flg2) { /* The odd nodes are evenly distributed across last nodes */ 17347c6ae99SBarry Smith x = (M + rank)/m; 1748865f1eaSKarl Rupp if (M/m == x) xs = rank*x; 1758865f1eaSKarl Rupp else xs = rank*(x-1) + (M+rank)%(x*m); 17647c6ae99SBarry Smith } else { /* The odd nodes are evenly distributed across the first k nodes */ 17747c6ae99SBarry Smith /* Regular PETSc Distribution */ 17847c6ae99SBarry Smith x = M/m + ((M % m) > rank); 1798865f1eaSKarl Rupp if (rank >= (M % m)) xs = (rank * (PetscInt)(M/m) + M % m); 1808865f1eaSKarl Rupp else xs = rank * (PetscInt)(M/m) + rank; 18147c6ae99SBarry Smith } 182fe16a2e9SBarry Smith ierr = MPI_Allgather(&xs,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr); 183fe16a2e9SBarry Smith for (i=0; i<m-1; i++) dd->lx[i] = dd->lx[i+1] - dd->lx[i]; 184fe16a2e9SBarry Smith dd->lx[m-1] = M - dd->lx[m-1]; 18547c6ae99SBarry Smith } else { 18647c6ae99SBarry Smith x = lx[rank]; 18747c6ae99SBarry Smith xs = 0; 1888865f1eaSKarl Rupp for (i=0; i<rank; i++) xs += lx[i]; 18947c6ae99SBarry Smith /* verify that data user provided is consistent */ 19047c6ae99SBarry Smith left = xs; 1918865f1eaSKarl Rupp for (i=rank; i<size; i++) left += lx[i]; 19247c6ae99SBarry Smith if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M %D %D",left,M); 19347c6ae99SBarry Smith } 19447c6ae99SBarry Smith 195bcea557cSEthan Coon /* 196bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 197bcea557cSEthan Coon the domain more than once 198bcea557cSEthan Coon */ 199bff4a2f0SMatthew G. Knepley if ((x < s) & ((M > 1) | (bx == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s); 200bcea557cSEthan Coon 20147c6ae99SBarry Smith xe = xs + x; 20247c6ae99SBarry Smith 20388661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 204d9c9ebe5SPeter Brune if (xs-sDist > 0) { 205d9c9ebe5SPeter Brune Xs = xs - sDist; 206d9c9ebe5SPeter Brune IXs = xs - sDist; 20788661749SPeter Brune } else { 2088865f1eaSKarl Rupp if (bx) Xs = xs - sDist; 2098865f1eaSKarl Rupp else Xs = 0; 21088661749SPeter Brune IXs = 0; 21188661749SPeter Brune } 21245b6f7e9SBarry Smith if (xe+sDist <= M) { 213d9c9ebe5SPeter Brune Xe = xe + sDist; 214d9c9ebe5SPeter Brune IXe = xe + sDist; 21588661749SPeter Brune } else { 2168865f1eaSKarl Rupp if (bx) Xe = xe + sDist; 21745b6f7e9SBarry Smith else Xe = M; 21845b6f7e9SBarry Smith IXe = M; 21947c6ae99SBarry Smith } 22047c6ae99SBarry Smith 221bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 222d9c9ebe5SPeter Brune Xs = xs - sDist; 223d9c9ebe5SPeter Brune Xe = xe + sDist; 224d9c9ebe5SPeter Brune IXs = xs - sDist; 225d9c9ebe5SPeter Brune IXe = xe + sDist; 226ce00eea3SSatish Balay } 227ce00eea3SSatish Balay 22847c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 22945b6f7e9SBarry Smith dd->Nlocal = dof*x; 230b1fb7eb7SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr); 23145b6f7e9SBarry Smith dd->nlocal = dof*(Xe-Xs); 232b1fb7eb7SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr); 23347c6ae99SBarry Smith 234bd1fc5aeSBarry Smith ierr = VecGetOwnershipRange(global,&start,NULL);CHKERRQ(ierr); 23547c6ae99SBarry Smith 23647c6ae99SBarry Smith /* Create Global to Local Vector Scatter Context */ 23747c6ae99SBarry Smith /* global to local must retrieve ghost points */ 23845b6f7e9SBarry Smith ierr = ISCreateStride(comm,dof*(IXe-IXs),dof*(IXs-Xs),1,&to);CHKERRQ(ierr); 23947c6ae99SBarry Smith 240854ce69bSBarry Smith ierr = PetscMalloc1(x+2*sDist,&idx);CHKERRQ(ierr); 2413bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)da,(x+2*(sDist))*sizeof(PetscInt));CHKERRQ(ierr); 24247c6ae99SBarry Smith 2438865f1eaSKarl Rupp for (i=0; i<IXs-Xs; i++) idx[i] = -1; /* prepend with -1s if needed for ghosted case*/ 24447c6ae99SBarry Smith 245ce00eea3SSatish Balay nn = IXs-Xs; 246bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) { /* Handle all cases with periodic first */ 247d9c9ebe5SPeter Brune for (i=0; i<sDist; i++) { /* Left ghost points */ 2488865f1eaSKarl Rupp if ((xs-sDist+i)>=0) idx[nn++] = xs-sDist+i; 24945b6f7e9SBarry Smith else idx[nn++] = M+(xs-sDist+i); 25047c6ae99SBarry Smith } 25147c6ae99SBarry Smith 2528865f1eaSKarl Rupp for (i=0; i<x; i++) idx [nn++] = xs + i; /* Non-ghost points */ 25347c6ae99SBarry Smith 254d9c9ebe5SPeter Brune for (i=0; i<sDist; i++) { /* Right ghost points */ 25545b6f7e9SBarry Smith if ((xe+i)<M) idx [nn++] = xe+i; 25645b6f7e9SBarry Smith else idx [nn++] = (xe+i) - M; 25747c6ae99SBarry Smith } 258bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { /* Handle all cases with periodic first */ 25945b6f7e9SBarry Smith for (i=0; i<(sDist); i++) { /* Left ghost points */ 26045b6f7e9SBarry Smith if ((xs-sDist+i)>=0) idx[nn++] = xs-sDist+i; 26145b6f7e9SBarry Smith else idx[nn++] = sDist - i; 26205900f5bSBarry Smith } 26305900f5bSBarry Smith 2648865f1eaSKarl Rupp for (i=0; i<x; i++) idx [nn++] = xs + i; /* Non-ghost points */ 26505900f5bSBarry Smith 26645b6f7e9SBarry Smith for (i=0; i<(sDist); i++) { /* Right ghost points */ 267662a7babSBarry Smith if ((xe+i)<M) idx[nn++] = xe+i; 268288e7d53SBarry Smith else idx[nn++] = M - (i + 2); 26905900f5bSBarry Smith } 27005900f5bSBarry Smith } else { /* Now do all cases with no periodicity */ 2718865f1eaSKarl Rupp if (0 <= xs-sDist) { 2728865f1eaSKarl Rupp for (i=0; i<sDist; i++) idx[nn++] = xs - sDist + i; 2738865f1eaSKarl Rupp } else { 2748865f1eaSKarl Rupp for (i=0; i<xs; i++) idx[nn++] = i; 2758865f1eaSKarl Rupp } 27647c6ae99SBarry Smith 2778865f1eaSKarl Rupp for (i=0; i<x; i++) idx [nn++] = xs + i; 27847c6ae99SBarry Smith 27945b6f7e9SBarry Smith if ((xe+sDist)<=M) { 2808865f1eaSKarl Rupp for (i=0; i<sDist; i++) idx[nn++]=xe+i; 2818865f1eaSKarl Rupp } else { 28245b6f7e9SBarry Smith for (i=xe; i<M; i++) idx[nn++]=i; 2838865f1eaSKarl Rupp } 28447c6ae99SBarry Smith } 28547c6ae99SBarry Smith 286b1fb7eb7SBarry Smith ierr = ISCreateBlock(comm,dof,nn-IXs+Xs,&idx[IXs-Xs],PETSC_USE_POINTER,&from);CHKERRQ(ierr); 2879448b7f1SJunchao Zhang ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 2883bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr); 289fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 290fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 291fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 292fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 29347c6ae99SBarry Smith 294d6e23781SBarry Smith dd->xs = dof*xs; dd->xe = dof*xe; dd->ys = 0; dd->ye = 1; dd->zs = 0; dd->ze = 1; 295d6e23781SBarry Smith dd->Xs = dof*Xs; dd->Xe = dof*Xe; dd->Ys = 0; dd->Ye = 1; dd->Zs = 0; dd->Ze = 1; 29647c6ae99SBarry Smith 29747c6ae99SBarry Smith dd->gtol = gtol; 298d6e23781SBarry Smith dd->base = dof*xs; 2999a42bb27SBarry Smith da->ops->view = DMView_DA_1d; 30047c6ae99SBarry Smith 30147c6ae99SBarry Smith /* 30247c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 30347c6ae99SBarry Smith of VecSetValuesLocal(). 30447c6ae99SBarry Smith */ 3058865f1eaSKarl Rupp for (i=0; i<Xe-IXe; i++) idx[nn++] = -1; /* pad with -1s if needed for ghosted case*/ 306ce00eea3SSatish Balay 30745b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr); 3083bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr); 30947c6ae99SBarry Smith 31047c6ae99SBarry Smith PetscFunctionReturn(0); 31147c6ae99SBarry Smith } 312d7bf68aeSBarry Smith 31347c6ae99SBarry Smith 31447c6ae99SBarry Smith /*@C 315aa219208SBarry Smith DMDACreate1d - Creates an object that will manage the communication of one-dimensional 31647c6ae99SBarry Smith regular array data that is distributed across some processors. 31747c6ae99SBarry Smith 318*d083f849SBarry Smith Collective 31947c6ae99SBarry Smith 32047c6ae99SBarry Smith Input Parameters: 32147c6ae99SBarry Smith + comm - MPI communicator 3221321219cSEthan Coon . bx - type of ghost cells at the boundary the array should have, if any. Use 323bff4a2f0SMatthew G. Knepley DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, or DM_BOUNDARY_PERIODIC. 324e43dc8daSBarry Smith . M - global dimension of the array (that is the number of grid points) 32547c6ae99SBarry Smith from the command line with -da_grid_x <M>) 32647c6ae99SBarry Smith . dof - number of degrees of freedom per node 32747c6ae99SBarry Smith . s - stencil width 32847c6ae99SBarry Smith - lx - array containing number of nodes in the X direction on each processor, 3290298fd71SBarry Smith or NULL. If non-null, must be of length as the number of processes in the MPI_Comm. 330e43dc8daSBarry Smith The sum of these entries must equal M 33147c6ae99SBarry Smith 33247c6ae99SBarry Smith Output Parameter: 33347c6ae99SBarry Smith . da - the resulting distributed array object 33447c6ae99SBarry Smith 33547c6ae99SBarry Smith Options Database Key: 336706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate1d() 337e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 338e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement factor 339e43dc8daSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it 34047c6ae99SBarry Smith 34147c6ae99SBarry Smith Level: beginner 34247c6ae99SBarry Smith 34347c6ae99SBarry Smith Notes: 344aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 345564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 346564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 34747c6ae99SBarry Smith 348897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 349897f7067SBarry Smith 350897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 351897f7067SBarry Smith but before DMSetUp(). 352897f7067SBarry Smith 353aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDASetRefinementFactor(), 35499f0b487SRichard Tran Mills DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDAGetRefinementFactor(), 355d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 35647c6ae99SBarry Smith 35747c6ae99SBarry Smith @*/ 358bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate1d(MPI_Comm comm, DMBoundaryType bx, PetscInt M, PetscInt dof, PetscInt s, const PetscInt lx[], DM *da) 35947c6ae99SBarry Smith { 36047c6ae99SBarry Smith PetscErrorCode ierr; 36147c6ae99SBarry Smith PetscMPIInt size; 36247c6ae99SBarry Smith 36347c6ae99SBarry Smith PetscFunctionBegin; 364aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 365c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*da, 1);CHKERRQ(ierr); 366aa219208SBarry Smith ierr = DMDASetSizes(*da, M, 1, 1);CHKERRQ(ierr); 36747c6ae99SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 368aa219208SBarry Smith ierr = DMDASetNumProcs(*da, size, PETSC_DECIDE, PETSC_DECIDE);CHKERRQ(ierr); 369bff4a2f0SMatthew G. Knepley ierr = DMDASetBoundaryType(*da, bx, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE);CHKERRQ(ierr); 370aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 371aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 3720298fd71SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, NULL, NULL);CHKERRQ(ierr); 37347c6ae99SBarry Smith PetscFunctionReturn(0); 37447c6ae99SBarry Smith } 375