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; 149a42bb27SBarry Smith PetscBool iascii,isdraw,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); 25251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 269a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 27251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr); 289a42bb27SBarry Smith #endif 2947c6ae99SBarry Smith if (iascii) { 3047c6ae99SBarry Smith PetscViewerFormat format; 3147c6ae99SBarry Smith 3247c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 33*8135c375SStefano Zampini if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL && format != PETSC_VIEWER_ASCII_GLVIS) { 34aa219208SBarry Smith DMDALocalInfo info; 35aa219208SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 361575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 3747c6ae99SBarry 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); 3847c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D\n",info.xs,info.xs+info.xm);CHKERRQ(ierr); 3947c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 401575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 41*8135c375SStefano Zampini } else if (format == PETSC_VIEWER_ASCII_GLVIS) { 42*8135c375SStefano Zampini ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr); 439a42bb27SBarry Smith } else { 449a42bb27SBarry Smith ierr = DMView_DA_VTK(da, viewer);CHKERRQ(ierr); 4547c6ae99SBarry Smith } 4647c6ae99SBarry Smith } else if (isdraw) { 4747c6ae99SBarry Smith PetscDraw draw; 4847c6ae99SBarry Smith double ymin = -1,ymax = 1,xmin = -1,xmax = dd->M,x; 4947c6ae99SBarry Smith PetscInt base; 5047c6ae99SBarry Smith char node[10]; 5147c6ae99SBarry Smith PetscBool isnull; 5247c6ae99SBarry Smith 5347c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 545b399a63SLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 555b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 565b399a63SLisandro Dalcin 575b399a63SLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 585b399a63SLisandro Dalcin ierr = PetscDrawClear(draw);CHKERRQ(ierr); 5947c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 6047c6ae99SBarry Smith 615b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 6247c6ae99SBarry Smith /* first processor draws all node lines */ 6347c6ae99SBarry Smith if (!rank) { 6447c6ae99SBarry Smith PetscInt xmin_tmp; 6547c6ae99SBarry Smith ymin = 0.0; ymax = 0.3; 6647c6ae99SBarry Smith for (xmin_tmp=0; xmin_tmp < dd->M; xmin_tmp++) { 6747c6ae99SBarry Smith ierr = PetscDrawLine(draw,(double)xmin_tmp,ymin,(double)xmin_tmp,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 6847c6ae99SBarry Smith } 6947c6ae99SBarry Smith xmin = 0.0; xmax = dd->M - 1; 7047c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7147c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7247c6ae99SBarry Smith } 735b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 745b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 7547c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 7647c6ae99SBarry Smith 775b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 7847c6ae99SBarry Smith /* draw my box */ 7947c6ae99SBarry Smith ymin = 0; ymax = 0.3; xmin = dd->xs / dd->w; xmax = (dd->xe / dd->w) - 1; 8047c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 8147c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 8247c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 8347c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 8447c6ae99SBarry Smith /* Put in index numbers */ 8547c6ae99SBarry Smith base = dd->base / dd->w; 8647c6ae99SBarry Smith for (x=xmin; x<=xmax; x++) { 87832b7cebSLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr); 8847c6ae99SBarry Smith ierr = PetscDrawString(draw,x,ymin,PETSC_DRAW_RED,node);CHKERRQ(ierr); 8947c6ae99SBarry Smith } 905b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 915b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 9247c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 93832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 949a42bb27SBarry Smith } else if (isbinary) { 959a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 969a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 979a42bb27SBarry Smith } else if (ismatlab) { 989a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 999a42bb27SBarry Smith #endif 10011aeaf0aSBarry Smith } 10147c6ae99SBarry Smith PetscFunctionReturn(0); 10247c6ae99SBarry Smith } 10347c6ae99SBarry Smith 10447c6ae99SBarry Smith 1057087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_1D(DM da) 10647c6ae99SBarry Smith { 10747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 10847c6ae99SBarry Smith const PetscInt M = dd->M; 10947c6ae99SBarry Smith const PetscInt dof = dd->w; 11047c6ae99SBarry Smith const PetscInt s = dd->s; 11145b6f7e9SBarry Smith const PetscInt sDist = s; /* stencil distance in points */ 11247c6ae99SBarry Smith const PetscInt *lx = dd->lx; 113bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 11447c6ae99SBarry Smith MPI_Comm comm; 11547c6ae99SBarry Smith Vec local, global; 116bd1fc5aeSBarry Smith VecScatter gtol; 11747c6ae99SBarry Smith IS to, from; 11847c6ae99SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE; 11947c6ae99SBarry Smith PetscMPIInt rank, size; 120bd1fc5aeSBarry Smith PetscInt i,*idx,nn,left,xs,xe,x,Xs,Xe,start,m,IXs,IXe; 12147c6ae99SBarry Smith PetscErrorCode ierr; 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith PetscFunctionBegin; 12447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 12547c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 12647c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 12747c6ae99SBarry Smith 1287d310018SBarry Smith dd->p = 1; 1297d310018SBarry Smith dd->n = 1; 13047c6ae99SBarry Smith dd->m = size; 13147c6ae99SBarry Smith m = dd->m; 13247c6ae99SBarry Smith 13347c6ae99SBarry Smith if (s > 0) { 13447c6ae99SBarry Smith /* if not communicating data then should be ok to have nothing on some processes */ 13547c6ae99SBarry Smith if (M < m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"More processes than data points! %D %D",m,M); 136a66092f1SBarry 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); 13747c6ae99SBarry Smith } 13847c6ae99SBarry Smith 13947c6ae99SBarry Smith /* 14047c6ae99SBarry Smith Determine locally owned region 14147c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 14247c6ae99SBarry Smith */ 14347c6ae99SBarry Smith if (!lx) { 144785e854fSJed Brown ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr); 145c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_partition_blockcomm",&flg1,NULL);CHKERRQ(ierr); 146c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_partition_nodes_at_end",&flg2,NULL);CHKERRQ(ierr); 14747c6ae99SBarry Smith if (flg1) { /* Block Comm type Distribution */ 14847c6ae99SBarry Smith xs = rank*M/m; 14947c6ae99SBarry Smith x = (rank + 1)*M/m - xs; 15047c6ae99SBarry Smith } else if (flg2) { /* The odd nodes are evenly distributed across last nodes */ 15147c6ae99SBarry Smith x = (M + rank)/m; 1528865f1eaSKarl Rupp if (M/m == x) xs = rank*x; 1538865f1eaSKarl Rupp else xs = rank*(x-1) + (M+rank)%(x*m); 15447c6ae99SBarry Smith } else { /* The odd nodes are evenly distributed across the first k nodes */ 15547c6ae99SBarry Smith /* Regular PETSc Distribution */ 15647c6ae99SBarry Smith x = M/m + ((M % m) > rank); 1578865f1eaSKarl Rupp if (rank >= (M % m)) xs = (rank * (PetscInt)(M/m) + M % m); 1588865f1eaSKarl Rupp else xs = rank * (PetscInt)(M/m) + rank; 15947c6ae99SBarry Smith } 160fe16a2e9SBarry Smith ierr = MPI_Allgather(&xs,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr); 161fe16a2e9SBarry Smith for (i=0; i<m-1; i++) dd->lx[i] = dd->lx[i+1] - dd->lx[i]; 162fe16a2e9SBarry Smith dd->lx[m-1] = M - dd->lx[m-1]; 16347c6ae99SBarry Smith } else { 16447c6ae99SBarry Smith x = lx[rank]; 16547c6ae99SBarry Smith xs = 0; 1668865f1eaSKarl Rupp for (i=0; i<rank; i++) xs += lx[i]; 16747c6ae99SBarry Smith /* verify that data user provided is consistent */ 16847c6ae99SBarry Smith left = xs; 1698865f1eaSKarl Rupp for (i=rank; i<size; i++) left += lx[i]; 17047c6ae99SBarry 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); 17147c6ae99SBarry Smith } 17247c6ae99SBarry Smith 173bcea557cSEthan Coon /* 174bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 175bcea557cSEthan Coon the domain more than once 176bcea557cSEthan Coon */ 177bff4a2f0SMatthew 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); 178bcea557cSEthan Coon 17947c6ae99SBarry Smith xe = xs + x; 18047c6ae99SBarry Smith 18188661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 182d9c9ebe5SPeter Brune if (xs-sDist > 0) { 183d9c9ebe5SPeter Brune Xs = xs - sDist; 184d9c9ebe5SPeter Brune IXs = xs - sDist; 18588661749SPeter Brune } else { 1868865f1eaSKarl Rupp if (bx) Xs = xs - sDist; 1878865f1eaSKarl Rupp else Xs = 0; 18888661749SPeter Brune IXs = 0; 18988661749SPeter Brune } 19045b6f7e9SBarry Smith if (xe+sDist <= M) { 191d9c9ebe5SPeter Brune Xe = xe + sDist; 192d9c9ebe5SPeter Brune IXe = xe + sDist; 19388661749SPeter Brune } else { 1948865f1eaSKarl Rupp if (bx) Xe = xe + sDist; 19545b6f7e9SBarry Smith else Xe = M; 19645b6f7e9SBarry Smith IXe = M; 19747c6ae99SBarry Smith } 19847c6ae99SBarry Smith 199bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 200d9c9ebe5SPeter Brune Xs = xs - sDist; 201d9c9ebe5SPeter Brune Xe = xe + sDist; 202d9c9ebe5SPeter Brune IXs = xs - sDist; 203d9c9ebe5SPeter Brune IXe = xe + sDist; 204ce00eea3SSatish Balay } 205ce00eea3SSatish Balay 20647c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 20745b6f7e9SBarry Smith dd->Nlocal = dof*x; 208b1fb7eb7SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr); 20945b6f7e9SBarry Smith dd->nlocal = dof*(Xe-Xs); 210b1fb7eb7SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr); 21147c6ae99SBarry Smith 212bd1fc5aeSBarry Smith ierr = VecGetOwnershipRange(global,&start,NULL);CHKERRQ(ierr); 21347c6ae99SBarry Smith 21447c6ae99SBarry Smith /* Create Global to Local Vector Scatter Context */ 21547c6ae99SBarry Smith /* global to local must retrieve ghost points */ 21645b6f7e9SBarry Smith ierr = ISCreateStride(comm,dof*(IXe-IXs),dof*(IXs-Xs),1,&to);CHKERRQ(ierr); 21747c6ae99SBarry Smith 218854ce69bSBarry Smith ierr = PetscMalloc1(x+2*sDist,&idx);CHKERRQ(ierr); 2193bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)da,(x+2*(sDist))*sizeof(PetscInt));CHKERRQ(ierr); 22047c6ae99SBarry Smith 2218865f1eaSKarl Rupp for (i=0; i<IXs-Xs; i++) idx[i] = -1; /* prepend with -1s if needed for ghosted case*/ 22247c6ae99SBarry Smith 223ce00eea3SSatish Balay nn = IXs-Xs; 224bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) { /* Handle all cases with periodic first */ 225d9c9ebe5SPeter Brune for (i=0; i<sDist; i++) { /* Left ghost points */ 2268865f1eaSKarl Rupp if ((xs-sDist+i)>=0) idx[nn++] = xs-sDist+i; 22745b6f7e9SBarry Smith else idx[nn++] = M+(xs-sDist+i); 22847c6ae99SBarry Smith } 22947c6ae99SBarry Smith 2308865f1eaSKarl Rupp for (i=0; i<x; i++) idx [nn++] = xs + i; /* Non-ghost points */ 23147c6ae99SBarry Smith 232d9c9ebe5SPeter Brune for (i=0; i<sDist; i++) { /* Right ghost points */ 23345b6f7e9SBarry Smith if ((xe+i)<M) idx [nn++] = xe+i; 23445b6f7e9SBarry Smith else idx [nn++] = (xe+i) - M; 23547c6ae99SBarry Smith } 236bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { /* Handle all cases with periodic first */ 23745b6f7e9SBarry Smith for (i=0; i<(sDist); i++) { /* Left ghost points */ 23845b6f7e9SBarry Smith if ((xs-sDist+i)>=0) idx[nn++] = xs-sDist+i; 23945b6f7e9SBarry Smith else idx[nn++] = sDist - i; 24005900f5bSBarry Smith } 24105900f5bSBarry Smith 2428865f1eaSKarl Rupp for (i=0; i<x; i++) idx [nn++] = xs + i; /* Non-ghost points */ 24305900f5bSBarry Smith 24445b6f7e9SBarry Smith for (i=0; i<(sDist); i++) { /* Right ghost points */ 245662a7babSBarry Smith if ((xe+i)<M) idx[nn++] = xe+i; 246288e7d53SBarry Smith else idx[nn++] = M - (i + 2); 24705900f5bSBarry Smith } 24805900f5bSBarry Smith } else { /* Now do all cases with no periodicity */ 2498865f1eaSKarl Rupp if (0 <= xs-sDist) { 2508865f1eaSKarl Rupp for (i=0; i<sDist; i++) idx[nn++] = xs - sDist + i; 2518865f1eaSKarl Rupp } else { 2528865f1eaSKarl Rupp for (i=0; i<xs; i++) idx[nn++] = i; 2538865f1eaSKarl Rupp } 25447c6ae99SBarry Smith 2558865f1eaSKarl Rupp for (i=0; i<x; i++) idx [nn++] = xs + i; 25647c6ae99SBarry Smith 25745b6f7e9SBarry Smith if ((xe+sDist)<=M) { 2588865f1eaSKarl Rupp for (i=0; i<sDist; i++) idx[nn++]=xe+i; 2598865f1eaSKarl Rupp } else { 26045b6f7e9SBarry Smith for (i=xe; i<M; i++) idx[nn++]=i; 2618865f1eaSKarl Rupp } 26247c6ae99SBarry Smith } 26347c6ae99SBarry Smith 264b1fb7eb7SBarry Smith ierr = ISCreateBlock(comm,dof,nn-IXs+Xs,&idx[IXs-Xs],PETSC_USE_POINTER,&from);CHKERRQ(ierr); 26547c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 2663bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr); 267fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 268fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 269fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 270fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 27147c6ae99SBarry Smith 272d6e23781SBarry Smith dd->xs = dof*xs; dd->xe = dof*xe; dd->ys = 0; dd->ye = 1; dd->zs = 0; dd->ze = 1; 273d6e23781SBarry Smith dd->Xs = dof*Xs; dd->Xe = dof*Xe; dd->Ys = 0; dd->Ye = 1; dd->Zs = 0; dd->Ze = 1; 27447c6ae99SBarry Smith 27547c6ae99SBarry Smith dd->gtol = gtol; 276d6e23781SBarry Smith dd->base = dof*xs; 2779a42bb27SBarry Smith da->ops->view = DMView_DA_1d; 27847c6ae99SBarry Smith 27947c6ae99SBarry Smith /* 28047c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 28147c6ae99SBarry Smith of VecSetValuesLocal(). 28247c6ae99SBarry Smith */ 2838865f1eaSKarl Rupp for (i=0; i<Xe-IXe; i++) idx[nn++] = -1; /* pad with -1s if needed for ghosted case*/ 284ce00eea3SSatish Balay 28545b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr); 2863bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr); 28747c6ae99SBarry Smith 28847c6ae99SBarry Smith PetscFunctionReturn(0); 28947c6ae99SBarry Smith } 290d7bf68aeSBarry Smith 29147c6ae99SBarry Smith 29247c6ae99SBarry Smith /*@C 293aa219208SBarry Smith DMDACreate1d - Creates an object that will manage the communication of one-dimensional 29447c6ae99SBarry Smith regular array data that is distributed across some processors. 29547c6ae99SBarry Smith 29647c6ae99SBarry Smith Collective on MPI_Comm 29747c6ae99SBarry Smith 29847c6ae99SBarry Smith Input Parameters: 29947c6ae99SBarry Smith + comm - MPI communicator 3001321219cSEthan Coon . bx - type of ghost cells at the boundary the array should have, if any. Use 301bff4a2f0SMatthew G. Knepley DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, or DM_BOUNDARY_PERIODIC. 302897f7067SBarry Smith . M - global dimension of the array 30347c6ae99SBarry Smith from the command line with -da_grid_x <M>) 30447c6ae99SBarry Smith . dof - number of degrees of freedom per node 30547c6ae99SBarry Smith . s - stencil width 30647c6ae99SBarry Smith - lx - array containing number of nodes in the X direction on each processor, 3070298fd71SBarry Smith or NULL. If non-null, must be of length as the number of processes in the MPI_Comm. 30847c6ae99SBarry Smith 30947c6ae99SBarry Smith Output Parameter: 31047c6ae99SBarry Smith . da - the resulting distributed array object 31147c6ae99SBarry Smith 31247c6ae99SBarry Smith Options Database Key: 313706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate1d() 31447c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction; can set if M < 0 315e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement factor 316e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it, if M < 0 31747c6ae99SBarry Smith 31847c6ae99SBarry Smith Level: beginner 31947c6ae99SBarry Smith 32047c6ae99SBarry Smith Notes: 321aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 322564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 323564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 32447c6ae99SBarry Smith 325897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 326897f7067SBarry Smith 327897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 328897f7067SBarry Smith but before DMSetUp(). 329897f7067SBarry Smith 33047c6ae99SBarry Smith .keywords: distributed array, create, one-dimensional 33147c6ae99SBarry Smith 332aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDASetRefinementFactor(), 33399f0b487SRichard Tran Mills DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDAGetRefinementFactor(), 334d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 33547c6ae99SBarry Smith 33647c6ae99SBarry Smith @*/ 337bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate1d(MPI_Comm comm, DMBoundaryType bx, PetscInt M, PetscInt dof, PetscInt s, const PetscInt lx[], DM *da) 33847c6ae99SBarry Smith { 33947c6ae99SBarry Smith PetscErrorCode ierr; 34047c6ae99SBarry Smith PetscMPIInt size; 34147c6ae99SBarry Smith 34247c6ae99SBarry Smith PetscFunctionBegin; 343aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 344c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*da, 1);CHKERRQ(ierr); 345aa219208SBarry Smith ierr = DMDASetSizes(*da, M, 1, 1);CHKERRQ(ierr); 34647c6ae99SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 347aa219208SBarry Smith ierr = DMDASetNumProcs(*da, size, PETSC_DECIDE, PETSC_DECIDE);CHKERRQ(ierr); 348bff4a2f0SMatthew G. Knepley ierr = DMDASetBoundaryType(*da, bx, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE);CHKERRQ(ierr); 349aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 350aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 3510298fd71SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, NULL, NULL);CHKERRQ(ierr); 35247c6ae99SBarry Smith PetscFunctionReturn(0); 35347c6ae99SBarry Smith } 354