xref: /petsc/src/dm/impls/da/da1.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
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   PetscMPIInt    rank;
13c9493c35SStefano Zampini   PetscBool      iascii,isdraw,isglvis,isbinary;
1447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
159a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
169a42bb27SBarry Smith   PetscBool ismatlab;
179a42bb27SBarry Smith #endif
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank));
2147c6ae99SBarry Smith 
229566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
239566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
249566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis));
259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
269a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab));
289a42bb27SBarry Smith #endif
2947c6ae99SBarry Smith   if (iascii) {
3047c6ae99SBarry Smith     PetscViewerFormat format;
3147c6ae99SBarry Smith 
329566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
3376a8abe0SBarry Smith     if (format == PETSC_VIEWER_LOAD_BALANCE) {
3476a8abe0SBarry Smith       PetscInt      i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal;
3576a8abe0SBarry Smith       DMDALocalInfo info;
3676a8abe0SBarry Smith       PetscMPIInt   size;
379566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da),&size));
389566063dSJacob Faibussowitsch       PetscCall(DMDAGetLocalInfo(da,&info));
3976a8abe0SBarry Smith       nzlocal = info.xm;
409566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(size,&nz));
419566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da)));
4276a8abe0SBarry Smith       for (i=0; i<(PetscInt)size; i++) {
4376a8abe0SBarry Smith         nmax = PetscMax(nmax,nz[i]);
4476a8abe0SBarry Smith         nmin = PetscMin(nmin,nz[i]);
4576a8abe0SBarry Smith         navg += nz[i];
4676a8abe0SBarry Smith       }
479566063dSJacob Faibussowitsch       PetscCall(PetscFree(nz));
4876a8abe0SBarry Smith       navg = navg/size;
4963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  Load Balance - Grid Points: Min %" PetscInt_FMT "  avg %" PetscInt_FMT "  max %" PetscInt_FMT "\n",nmin,navg,nmax));
5076a8abe0SBarry Smith       PetscFunctionReturn(0);
5176a8abe0SBarry Smith     }
528ec8862eSJed Brown     if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) {
53aa219208SBarry Smith       DMDALocalInfo info;
549566063dSJacob Faibussowitsch       PetscCall(DMDAGetLocalInfo(da,&info));
559566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
5663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %" PetscInt_FMT " m %" PetscInt_FMT " w %" PetscInt_FMT " s %" PetscInt_FMT "\n",rank,dd->M,dd->m,dd->w,dd->s));
5763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %" PetscInt_FMT " %" PetscInt_FMT "\n",info.xs,info.xs+info.xm));
589566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
599566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
60*1baa6e33SBarry Smith     } else if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMView_DA_GLVis(da,viewer));
61*1baa6e33SBarry Smith     else PetscCall(DMView_DA_VTK(da, viewer));
6247c6ae99SBarry Smith   } else if (isdraw) {
6347c6ae99SBarry Smith     PetscDraw      draw;
6447c6ae99SBarry Smith     double         ymin = -1,ymax = 1,xmin = -1,xmax = dd->M,x;
6547c6ae99SBarry Smith     PetscInt       base;
6647c6ae99SBarry Smith     char           node[10];
6747c6ae99SBarry Smith     PetscBool      isnull;
6847c6ae99SBarry Smith 
699566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
709566063dSJacob Faibussowitsch     PetscCall(PetscDrawIsNull(draw,&isnull));
715b399a63SLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
725b399a63SLisandro Dalcin 
739566063dSJacob Faibussowitsch     PetscCall(PetscDrawCheckResizedWindow(draw));
749566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
759566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax));
7647c6ae99SBarry Smith 
77d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
7847c6ae99SBarry Smith     /* first processor draws all node lines */
79dd400576SPatrick Sanan     if (rank == 0) {
8047c6ae99SBarry Smith       PetscInt xmin_tmp;
8147c6ae99SBarry Smith       ymin = 0.0; ymax = 0.3;
8247c6ae99SBarry Smith       for (xmin_tmp=0; xmin_tmp < dd->M; xmin_tmp++) {
839566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,(double)xmin_tmp,ymin,(double)xmin_tmp,ymax,PETSC_DRAW_BLACK));
8447c6ae99SBarry Smith       }
8547c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
869566063dSJacob Faibussowitsch       PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK));
879566063dSJacob Faibussowitsch       PetscCall(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_BLACK));
8847c6ae99SBarry Smith     }
89d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
909566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
919566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
9247c6ae99SBarry Smith 
93d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
9447c6ae99SBarry Smith     /* draw my box */
9547c6ae99SBarry Smith     ymin = 0; ymax = 0.3; xmin = dd->xs / dd->w; xmax = (dd->xe / dd->w)  - 1;
969566063dSJacob Faibussowitsch     PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED));
979566063dSJacob Faibussowitsch     PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED));
989566063dSJacob Faibussowitsch     PetscCall(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED));
999566063dSJacob Faibussowitsch     PetscCall(PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED));
10047c6ae99SBarry Smith     /* Put in index numbers */
10147c6ae99SBarry Smith     base = dd->base / dd->w;
10247c6ae99SBarry Smith     for (x=xmin; x<=xmax; x++) {
1039566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)base++));
1049566063dSJacob Faibussowitsch       PetscCall(PetscDrawString(draw,x,ymin,PETSC_DRAW_RED,node));
10547c6ae99SBarry Smith     }
106d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1079566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
1089566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
1099566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(draw));
110c9493c35SStefano Zampini   } else if (isglvis) {
1119566063dSJacob Faibussowitsch     PetscCall(DMView_DA_GLVis(da,viewer));
1129a42bb27SBarry Smith   } else if (isbinary) {
1139566063dSJacob Faibussowitsch     PetscCall(DMView_DA_Binary(da,viewer));
1149a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1159a42bb27SBarry Smith   } else if (ismatlab) {
1169566063dSJacob Faibussowitsch     PetscCall(DMView_DA_Matlab(da,viewer));
1179a42bb27SBarry Smith #endif
11811aeaf0aSBarry Smith   }
11947c6ae99SBarry Smith   PetscFunctionReturn(0);
12047c6ae99SBarry Smith }
12147c6ae99SBarry Smith 
1227087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_1D(DM da)
12347c6ae99SBarry Smith {
12447c6ae99SBarry Smith   DM_DA            *dd   = (DM_DA*)da->data;
12547c6ae99SBarry Smith   const PetscInt   M     = dd->M;
12647c6ae99SBarry Smith   const PetscInt   dof   = dd->w;
12747c6ae99SBarry Smith   const PetscInt   s     = dd->s;
12845b6f7e9SBarry Smith   const PetscInt   sDist = s;  /* stencil distance in points */
12947c6ae99SBarry Smith   const PetscInt   *lx   = dd->lx;
130bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx    = dd->bx;
13147c6ae99SBarry Smith   MPI_Comm         comm;
13247c6ae99SBarry Smith   Vec              local, global;
133bd1fc5aeSBarry Smith   VecScatter       gtol;
13447c6ae99SBarry Smith   IS               to, from;
13547c6ae99SBarry Smith   PetscBool        flg1 = PETSC_FALSE, flg2 = PETSC_FALSE;
13647c6ae99SBarry Smith   PetscMPIInt      rank, size;
137bd1fc5aeSBarry Smith   PetscInt         i,*idx,nn,left,xs,xe,x,Xs,Xe,start,m,IXs,IXe;
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) da, &comm));
1419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
1429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
14347c6ae99SBarry Smith 
1447d310018SBarry Smith   dd->p = 1;
1457d310018SBarry Smith   dd->n = 1;
14647c6ae99SBarry Smith   dd->m = size;
14747c6ae99SBarry Smith   m     = dd->m;
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith   if (s > 0) {
15047c6ae99SBarry Smith     /* if not communicating data then should be ok to have nothing on some processes */
15163a3b9bcSJacob Faibussowitsch     PetscCheck(M >= m,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"More processes than data points! %" PetscInt_FMT " %" PetscInt_FMT,m,M);
15263a3b9bcSJacob Faibussowitsch     PetscCheck((M-1) >= s || size <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array is too small for stencil! %" PetscInt_FMT " %" PetscInt_FMT,M-1,s);
15347c6ae99SBarry Smith   }
15447c6ae99SBarry Smith 
15547c6ae99SBarry Smith   /*
15647c6ae99SBarry Smith      Determine locally owned region
15747c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
15847c6ae99SBarry Smith   */
15947c6ae99SBarry Smith   if (!lx) {
1609566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &dd->lx));
1619566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_partition_blockcomm",&flg1,NULL));
1629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_partition_nodes_at_end",&flg2,NULL));
16347c6ae99SBarry Smith     if (flg1) {      /* Block Comm type Distribution */
16447c6ae99SBarry Smith       xs = rank*M/m;
16547c6ae99SBarry Smith       x  = (rank + 1)*M/m - xs;
16647c6ae99SBarry Smith     } else if (flg2) { /* The odd nodes are evenly distributed across last nodes */
16747c6ae99SBarry Smith       x = (M + rank)/m;
1688865f1eaSKarl Rupp       if (M/m == x) xs = rank*x;
1698865f1eaSKarl Rupp       else          xs = rank*(x-1) + (M+rank)%(x*m);
17047c6ae99SBarry Smith     } else { /* The odd nodes are evenly distributed across the first k nodes */
17147c6ae99SBarry Smith       /* Regular PETSc Distribution */
17247c6ae99SBarry Smith       x = M/m + ((M % m) > rank);
1738865f1eaSKarl Rupp       if (rank >= (M % m)) xs = (rank * (PetscInt)(M/m) + M % m);
1748865f1eaSKarl Rupp       else                 xs = rank * (PetscInt)(M/m) + rank;
17547c6ae99SBarry Smith     }
1769566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgather(&xs,1,MPIU_INT,dd->lx,1,MPIU_INT,comm));
177fe16a2e9SBarry Smith     for (i=0; i<m-1; i++) dd->lx[i] = dd->lx[i+1] - dd->lx[i];
178fe16a2e9SBarry Smith     dd->lx[m-1] = M - dd->lx[m-1];
17947c6ae99SBarry Smith   } else {
18047c6ae99SBarry Smith     x  = lx[rank];
18147c6ae99SBarry Smith     xs = 0;
1828865f1eaSKarl Rupp     for (i=0; i<rank; i++) xs += lx[i];
18347c6ae99SBarry Smith     /* verify that data user provided is consistent */
18447c6ae99SBarry Smith     left = xs;
1858865f1eaSKarl Rupp     for (i=rank; i<size; i++) left += lx[i];
18663a3b9bcSJacob Faibussowitsch     PetscCheck(left == M,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M %" PetscInt_FMT " %" PetscInt_FMT,left,M);
18747c6ae99SBarry Smith   }
18847c6ae99SBarry Smith 
189bcea557cSEthan Coon   /*
190bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
191bcea557cSEthan Coon    the domain more than once
192bcea557cSEthan Coon   */
19363a3b9bcSJacob Faibussowitsch   PetscCheck((x >= s) || ((M <= 1) && (bx != DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT,x,s);
194bcea557cSEthan Coon 
19547c6ae99SBarry Smith   xe  = xs + x;
19647c6ae99SBarry Smith 
19788661749SPeter Brune   /* determine ghost region (Xs) and region scattered into (IXs)  */
198d9c9ebe5SPeter Brune   if (xs-sDist > 0) {
199d9c9ebe5SPeter Brune     Xs  = xs - sDist;
200d9c9ebe5SPeter Brune     IXs = xs - sDist;
20188661749SPeter Brune   } else {
2028865f1eaSKarl Rupp     if (bx) Xs = xs - sDist;
2038865f1eaSKarl Rupp     else Xs = 0;
20488661749SPeter Brune     IXs = 0;
20588661749SPeter Brune   }
20645b6f7e9SBarry Smith   if (xe+sDist <= M) {
207d9c9ebe5SPeter Brune     Xe  = xe + sDist;
208d9c9ebe5SPeter Brune     IXe = xe + sDist;
20988661749SPeter Brune   } else {
2108865f1eaSKarl Rupp     if (bx) Xe = xe + sDist;
21145b6f7e9SBarry Smith     else Xe = M;
21245b6f7e9SBarry Smith     IXe = M;
21347c6ae99SBarry Smith   }
21447c6ae99SBarry Smith 
215bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
216d9c9ebe5SPeter Brune     Xs  = xs - sDist;
217d9c9ebe5SPeter Brune     Xe  = xe + sDist;
218d9c9ebe5SPeter Brune     IXs = xs - sDist;
219d9c9ebe5SPeter Brune     IXe = xe + sDist;
220ce00eea3SSatish Balay   }
221ce00eea3SSatish Balay 
22247c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
22345b6f7e9SBarry Smith   dd->Nlocal = dof*x;
2249566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global));
22545b6f7e9SBarry Smith   dd->nlocal = dof*(Xe-Xs);
2269566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local));
22747c6ae99SBarry Smith 
2289566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(global,&start,NULL));
22947c6ae99SBarry Smith 
23047c6ae99SBarry Smith   /* Create Global to Local Vector Scatter Context */
23147c6ae99SBarry Smith   /* global to local must retrieve ghost points */
2329566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(comm,dof*(IXe-IXs),dof*(IXs-Xs),1,&to));
23347c6ae99SBarry Smith 
2349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(x+2*sDist,&idx));
2359566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)da,(x+2*(sDist))*sizeof(PetscInt)));
23647c6ae99SBarry Smith 
2378865f1eaSKarl Rupp   for (i=0; i<IXs-Xs; i++) idx[i] = -1; /* prepend with -1s if needed for ghosted case*/
23847c6ae99SBarry Smith 
239ce00eea3SSatish Balay   nn = IXs-Xs;
240bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) { /* Handle all cases with periodic first */
241d9c9ebe5SPeter Brune     for (i=0; i<sDist; i++) {  /* Left ghost points */
2428865f1eaSKarl Rupp       if ((xs-sDist+i)>=0) idx[nn++] = xs-sDist+i;
24345b6f7e9SBarry Smith       else                 idx[nn++] = M+(xs-sDist+i);
24447c6ae99SBarry Smith     }
24547c6ae99SBarry Smith 
2468865f1eaSKarl Rupp     for (i=0; i<x; i++) idx [nn++] = xs + i;  /* Non-ghost points */
24747c6ae99SBarry Smith 
248d9c9ebe5SPeter Brune     for (i=0; i<sDist; i++) { /* Right ghost points */
24945b6f7e9SBarry Smith       if ((xe+i)<M) idx [nn++] =  xe+i;
25045b6f7e9SBarry Smith       else          idx [nn++] = (xe+i) - M;
25147c6ae99SBarry Smith     }
252bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_MIRROR) { /* Handle all cases with periodic first */
25345b6f7e9SBarry Smith     for (i=0; i<(sDist); i++) {  /* Left ghost points */
25445b6f7e9SBarry Smith       if ((xs-sDist+i)>=0) idx[nn++] = xs-sDist+i;
25545b6f7e9SBarry Smith       else                 idx[nn++] = sDist - i;
25605900f5bSBarry Smith     }
25705900f5bSBarry Smith 
2588865f1eaSKarl Rupp     for (i=0; i<x; i++) idx [nn++] = xs + i;  /* Non-ghost points */
25905900f5bSBarry Smith 
26045b6f7e9SBarry Smith     for (i=0; i<(sDist); i++) { /* Right ghost points */
261662a7babSBarry Smith       if ((xe+i)<M) idx[nn++] =  xe+i;
262288e7d53SBarry Smith       else          idx[nn++] = M - (i + 2);
26305900f5bSBarry Smith     }
26405900f5bSBarry Smith   } else {      /* Now do all cases with no periodicity */
2658865f1eaSKarl Rupp     if (0 <= xs-sDist) {
2668865f1eaSKarl Rupp       for (i=0; i<sDist; i++) idx[nn++] = xs - sDist + i;
2678865f1eaSKarl Rupp     } else {
2688865f1eaSKarl Rupp       for (i=0; i<xs; i++) idx[nn++] = i;
2698865f1eaSKarl Rupp     }
27047c6ae99SBarry Smith 
2718865f1eaSKarl Rupp     for (i=0; i<x; i++) idx [nn++] = xs + i;
27247c6ae99SBarry Smith 
27345b6f7e9SBarry Smith     if ((xe+sDist)<=M) {
2748865f1eaSKarl Rupp       for (i=0; i<sDist; i++) idx[nn++]=xe+i;
2758865f1eaSKarl Rupp     } else {
27645b6f7e9SBarry Smith       for (i=xe; i<M; i++) idx[nn++]=i;
2778865f1eaSKarl Rupp     }
27847c6ae99SBarry Smith   }
27947c6ae99SBarry Smith 
2809566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(comm,dof,nn-IXs+Xs,&idx[IXs-Xs],PETSC_USE_POINTER,&from));
2819566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global,from,local,to,&gtol));
2829566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)gtol));
2839566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
2849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
2859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
2869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
28747c6ae99SBarry Smith 
288d6e23781SBarry Smith   dd->xs = dof*xs; dd->xe = dof*xe; dd->ys = 0; dd->ye = 1; dd->zs = 0; dd->ze = 1;
289d6e23781SBarry Smith   dd->Xs = dof*Xs; dd->Xe = dof*Xe; dd->Ys = 0; dd->Ye = 1; dd->Zs = 0; dd->Ze = 1;
29047c6ae99SBarry Smith 
29147c6ae99SBarry Smith   dd->gtol      = gtol;
292d6e23781SBarry Smith   dd->base      = dof*xs;
2939a42bb27SBarry Smith   da->ops->view = DMView_DA_1d;
29447c6ae99SBarry Smith 
29547c6ae99SBarry Smith   /*
29647c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
29747c6ae99SBarry Smith      of VecSetValuesLocal().
29847c6ae99SBarry Smith   */
2998865f1eaSKarl Rupp   for (i=0; i<Xe-IXe; i++) idx[nn++] = -1; /* pad with -1s if needed for ghosted case*/
300ce00eea3SSatish Balay 
3019566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap));
3029566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap));
30347c6ae99SBarry Smith 
30447c6ae99SBarry Smith   PetscFunctionReturn(0);
30547c6ae99SBarry Smith }
306d7bf68aeSBarry Smith 
30747c6ae99SBarry Smith /*@C
308aa219208SBarry Smith    DMDACreate1d - Creates an object that will manage the communication of  one-dimensional
30947c6ae99SBarry Smith    regular array data that is distributed across some processors.
31047c6ae99SBarry Smith 
311d083f849SBarry Smith    Collective
31247c6ae99SBarry Smith 
31347c6ae99SBarry Smith    Input Parameters:
31447c6ae99SBarry Smith +  comm - MPI communicator
3151321219cSEthan Coon .  bx - type of ghost cells at the boundary the array should have, if any. Use
316bff4a2f0SMatthew G. Knepley           DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, or DM_BOUNDARY_PERIODIC.
317e43dc8daSBarry Smith .  M - global dimension of the array (that is the number of grid points)
31847c6ae99SBarry Smith             from the command line with -da_grid_x <M>)
31947c6ae99SBarry Smith .  dof - number of degrees of freedom per node
32047c6ae99SBarry Smith .  s - stencil width
32147c6ae99SBarry Smith -  lx - array containing number of nodes in the X direction on each processor,
3220298fd71SBarry Smith         or NULL. If non-null, must be of length as the number of processes in the MPI_Comm.
323e43dc8daSBarry Smith         The sum of these entries must equal M
32447c6ae99SBarry Smith 
32547c6ae99SBarry Smith    Output Parameter:
32647c6ae99SBarry Smith .  da - the resulting distributed array object
32747c6ae99SBarry Smith 
32847c6ae99SBarry Smith    Options Database Key:
329706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate1d()
330e43dc8daSBarry Smith .  -da_grid_x <nx> - number of grid points in x direction
331e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement factor
332e43dc8daSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating it
33347c6ae99SBarry Smith 
33447c6ae99SBarry Smith    Level: beginner
33547c6ae99SBarry Smith 
33647c6ae99SBarry Smith    Notes:
337aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
338564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
339564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
34047c6ae99SBarry Smith 
341897f7067SBarry Smith    You must call DMSetUp() after this call before using this DM.
342897f7067SBarry Smith 
343897f7067SBarry Smith    If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call
344897f7067SBarry Smith    but before DMSetUp().
345897f7067SBarry Smith 
346db781477SPatrick Sanan .seealso: `DMDestroy()`, `DMView()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMGlobalToLocalBegin()`, `DMDASetRefinementFactor()`,
347db781477SPatrick Sanan           `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDAGetRefinementFactor()`,
348db781477SPatrick Sanan           `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`,
349db781477SPatrick Sanan           `DMStagCreate1d()`
35047c6ae99SBarry Smith 
35147c6ae99SBarry Smith @*/
352bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate1d(MPI_Comm comm, DMBoundaryType bx, PetscInt M, PetscInt dof, PetscInt s, const PetscInt lx[], DM *da)
35347c6ae99SBarry Smith {
35447c6ae99SBarry Smith   PetscMPIInt    size;
35547c6ae99SBarry Smith 
35647c6ae99SBarry Smith   PetscFunctionBegin;
3579566063dSJacob Faibussowitsch   PetscCall(DMDACreate(comm, da));
3589566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*da, 1));
3599566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(*da, M, 1, 1));
3609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3619566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(*da, size, PETSC_DECIDE, PETSC_DECIDE));
3629566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(*da, bx, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE));
3639566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(*da, dof));
3649566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(*da, s));
3659566063dSJacob Faibussowitsch   PetscCall(DMDASetOwnershipRanges(*da, lx, NULL, NULL));
36647c6ae99SBarry Smith   PetscFunctionReturn(0);
36747c6ae99SBarry Smith }
368