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> 10d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_DA_1d(DM da, PetscViewer viewer) 11d71ae5a4SJacob Faibussowitsch { 1247c6ae99SBarry Smith PetscMPIInt rank; 13c9493c35SStefano Zampini PetscBool iascii, isdraw, isglvis, isbinary; 1447c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 15d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 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)); 26d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 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)); 601baa6e33SBarry Smith } else if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMView_DA_GLVis(da, viewer)); 611baa6e33SBarry 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; 819371c9d4SSatish Balay ymin = 0.0; 829371c9d4SSatish Balay ymax = 0.3; 8348a46eb9SPierre Jolivet for (xmin_tmp = 0; xmin_tmp < dd->M; xmin_tmp++) PetscCall(PetscDrawLine(draw, (double)xmin_tmp, ymin, (double)xmin_tmp, ymax, PETSC_DRAW_BLACK)); 849371c9d4SSatish Balay xmin = 0.0; 859371c9d4SSatish Balay 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 */ 959371c9d4SSatish Balay ymin = 0; 969371c9d4SSatish Balay ymax = 0.3; 979371c9d4SSatish Balay xmin = dd->xs / dd->w; 989371c9d4SSatish Balay xmax = (dd->xe / dd->w) - 1; 999566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ymin, PETSC_DRAW_RED)); 1009566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymin, xmin, ymax, PETSC_DRAW_RED)); 1019566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymax, xmax, ymax, PETSC_DRAW_RED)); 1029566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmax, ymin, xmax, ymax, PETSC_DRAW_RED)); 10347c6ae99SBarry Smith /* Put in index numbers */ 10447c6ae99SBarry Smith base = dd->base / dd->w; 10547c6ae99SBarry Smith for (x = xmin; x <= xmax; x++) { 1069566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node, sizeof(node), "%d", (int)base++)); 1079566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, x, ymin, PETSC_DRAW_RED, node)); 10847c6ae99SBarry Smith } 109d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1109566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1119566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1129566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 113c9493c35SStefano Zampini } else if (isglvis) { 1149566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da, viewer)); 1159a42bb27SBarry Smith } else if (isbinary) { 1169566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da, viewer)); 117d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 1189a42bb27SBarry Smith } else if (ismatlab) { 1199566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da, viewer)); 1209a42bb27SBarry Smith #endif 12111aeaf0aSBarry Smith } 12247c6ae99SBarry Smith PetscFunctionReturn(0); 12347c6ae99SBarry Smith } 12447c6ae99SBarry Smith 125d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUp_DA_1D(DM da) 126d71ae5a4SJacob Faibussowitsch { 12747c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 12847c6ae99SBarry Smith const PetscInt M = dd->M; 12947c6ae99SBarry Smith const PetscInt dof = dd->w; 13047c6ae99SBarry Smith const PetscInt s = dd->s; 13145b6f7e9SBarry Smith const PetscInt sDist = s; /* stencil distance in points */ 13247c6ae99SBarry Smith const PetscInt *lx = dd->lx; 133bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 13447c6ae99SBarry Smith MPI_Comm comm; 13547c6ae99SBarry Smith Vec local, global; 136bd1fc5aeSBarry Smith VecScatter gtol; 13747c6ae99SBarry Smith IS to, from; 13847c6ae99SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE; 13947c6ae99SBarry Smith PetscMPIInt rank, size; 140bd1fc5aeSBarry Smith PetscInt i, *idx, nn, left, xs, xe, x, Xs, Xe, start, m, IXs, IXe; 14147c6ae99SBarry Smith 14247c6ae99SBarry Smith PetscFunctionBegin; 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 1449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 14647c6ae99SBarry Smith 1477d310018SBarry Smith dd->p = 1; 1487d310018SBarry Smith dd->n = 1; 14947c6ae99SBarry Smith dd->m = size; 15047c6ae99SBarry Smith m = dd->m; 15147c6ae99SBarry Smith 15247c6ae99SBarry Smith if (s > 0) { 15347c6ae99SBarry Smith /* if not communicating data then should be ok to have nothing on some processes */ 15463a3b9bcSJacob Faibussowitsch PetscCheck(M >= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "More processes than data points! %" PetscInt_FMT " %" PetscInt_FMT, m, M); 15563a3b9bcSJacob 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); 15647c6ae99SBarry Smith } 15747c6ae99SBarry Smith 15847c6ae99SBarry Smith /* 15947c6ae99SBarry Smith Determine locally owned region 16047c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 16147c6ae99SBarry Smith */ 16247c6ae99SBarry Smith if (!lx) { 1639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 1649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_partition_blockcomm", &flg1, NULL)); 1659566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_partition_nodes_at_end", &flg2, NULL)); 16647c6ae99SBarry Smith if (flg1) { /* Block Comm type Distribution */ 16747c6ae99SBarry Smith xs = rank * M / m; 16847c6ae99SBarry Smith x = (rank + 1) * M / m - xs; 16947c6ae99SBarry Smith } else if (flg2) { /* The odd nodes are evenly distributed across last nodes */ 17047c6ae99SBarry Smith x = (M + rank) / m; 1718865f1eaSKarl Rupp if (M / m == x) xs = rank * x; 1728865f1eaSKarl Rupp else xs = rank * (x - 1) + (M + rank) % (x * m); 17347c6ae99SBarry Smith } else { /* The odd nodes are evenly distributed across the first k nodes */ 17447c6ae99SBarry Smith /* Regular PETSc Distribution */ 17547c6ae99SBarry Smith x = M / m + ((M % m) > rank); 1768865f1eaSKarl Rupp if (rank >= (M % m)) xs = (rank * (PetscInt)(M / m) + M % m); 1778865f1eaSKarl Rupp else xs = rank * (PetscInt)(M / m) + rank; 17847c6ae99SBarry Smith } 1799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&xs, 1, MPIU_INT, dd->lx, 1, MPIU_INT, comm)); 180fe16a2e9SBarry Smith for (i = 0; i < m - 1; i++) dd->lx[i] = dd->lx[i + 1] - dd->lx[i]; 181fe16a2e9SBarry Smith dd->lx[m - 1] = M - dd->lx[m - 1]; 18247c6ae99SBarry Smith } else { 18347c6ae99SBarry Smith x = lx[rank]; 18447c6ae99SBarry Smith xs = 0; 1858865f1eaSKarl Rupp for (i = 0; i < rank; i++) xs += lx[i]; 18647c6ae99SBarry Smith /* verify that data user provided is consistent */ 18747c6ae99SBarry Smith left = xs; 1888865f1eaSKarl Rupp for (i = rank; i < size; i++) left += lx[i]; 18963a3b9bcSJacob 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); 19047c6ae99SBarry Smith } 19147c6ae99SBarry Smith 192bcea557cSEthan Coon /* 193bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 194bcea557cSEthan Coon the domain more than once 195bcea557cSEthan Coon */ 19663a3b9bcSJacob 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); 197bcea557cSEthan Coon 19847c6ae99SBarry Smith xe = xs + x; 19947c6ae99SBarry Smith 20088661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 201d9c9ebe5SPeter Brune if (xs - sDist > 0) { 202d9c9ebe5SPeter Brune Xs = xs - sDist; 203d9c9ebe5SPeter Brune IXs = xs - sDist; 20488661749SPeter Brune } else { 2058865f1eaSKarl Rupp if (bx) Xs = xs - sDist; 2068865f1eaSKarl Rupp else Xs = 0; 20788661749SPeter Brune IXs = 0; 20888661749SPeter Brune } 20945b6f7e9SBarry Smith if (xe + sDist <= M) { 210d9c9ebe5SPeter Brune Xe = xe + sDist; 211d9c9ebe5SPeter Brune IXe = xe + sDist; 21288661749SPeter Brune } else { 2138865f1eaSKarl Rupp if (bx) Xe = xe + sDist; 21445b6f7e9SBarry Smith else Xe = M; 21545b6f7e9SBarry Smith IXe = M; 21647c6ae99SBarry Smith } 21747c6ae99SBarry Smith 218bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 219d9c9ebe5SPeter Brune Xs = xs - sDist; 220d9c9ebe5SPeter Brune Xe = xe + sDist; 221d9c9ebe5SPeter Brune IXs = xs - sDist; 222d9c9ebe5SPeter Brune IXe = xe + sDist; 223ce00eea3SSatish Balay } 224ce00eea3SSatish Balay 22547c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 22645b6f7e9SBarry Smith dd->Nlocal = dof * x; 2279566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, dof, dd->Nlocal, PETSC_DECIDE, NULL, &global)); 22845b6f7e9SBarry Smith dd->nlocal = dof * (Xe - Xs); 2299566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dof, dd->nlocal, NULL, &local)); 23047c6ae99SBarry Smith 2319566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(global, &start, NULL)); 23247c6ae99SBarry Smith 23347c6ae99SBarry Smith /* Create Global to Local Vector Scatter Context */ 23447c6ae99SBarry Smith /* global to local must retrieve ghost points */ 2359566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, dof * (IXe - IXs), dof * (IXs - Xs), 1, &to)); 23647c6ae99SBarry Smith 2379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(x + 2 * sDist, &idx)); 23847c6ae99SBarry Smith 2398865f1eaSKarl Rupp for (i = 0; i < IXs - Xs; i++) idx[i] = -1; /* prepend with -1s if needed for ghosted case*/ 24047c6ae99SBarry Smith 241ce00eea3SSatish Balay nn = IXs - Xs; 242bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) { /* Handle all cases with periodic first */ 243d9c9ebe5SPeter Brune for (i = 0; i < sDist; i++) { /* Left ghost points */ 2448865f1eaSKarl Rupp if ((xs - sDist + i) >= 0) idx[nn++] = xs - sDist + i; 24545b6f7e9SBarry Smith else idx[nn++] = M + (xs - sDist + i); 24647c6ae99SBarry Smith } 24747c6ae99SBarry Smith 2488865f1eaSKarl Rupp for (i = 0; i < x; i++) idx[nn++] = xs + i; /* Non-ghost points */ 24947c6ae99SBarry Smith 250d9c9ebe5SPeter Brune for (i = 0; i < sDist; i++) { /* Right ghost points */ 25145b6f7e9SBarry Smith if ((xe + i) < M) idx[nn++] = xe + i; 25245b6f7e9SBarry Smith else idx[nn++] = (xe + i) - M; 25347c6ae99SBarry Smith } 254bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { /* Handle all cases with periodic first */ 25545b6f7e9SBarry Smith for (i = 0; i < (sDist); i++) { /* Left ghost points */ 25645b6f7e9SBarry Smith if ((xs - sDist + i) >= 0) idx[nn++] = xs - sDist + i; 25745b6f7e9SBarry Smith else idx[nn++] = sDist - i; 25805900f5bSBarry Smith } 25905900f5bSBarry Smith 2608865f1eaSKarl Rupp for (i = 0; i < x; i++) idx[nn++] = xs + i; /* Non-ghost points */ 26105900f5bSBarry Smith 26245b6f7e9SBarry Smith for (i = 0; i < (sDist); i++) { /* Right ghost points */ 263662a7babSBarry Smith if ((xe + i) < M) idx[nn++] = xe + i; 264288e7d53SBarry Smith else idx[nn++] = M - (i + 2); 26505900f5bSBarry Smith } 26605900f5bSBarry Smith } else { /* Now do all cases with no periodicity */ 2678865f1eaSKarl Rupp if (0 <= xs - sDist) { 2688865f1eaSKarl Rupp for (i = 0; i < sDist; i++) idx[nn++] = xs - sDist + i; 2698865f1eaSKarl Rupp } else { 2708865f1eaSKarl Rupp for (i = 0; i < xs; i++) idx[nn++] = i; 2718865f1eaSKarl Rupp } 27247c6ae99SBarry Smith 2738865f1eaSKarl Rupp for (i = 0; i < x; i++) idx[nn++] = xs + i; 27447c6ae99SBarry Smith 27545b6f7e9SBarry Smith if ((xe + sDist) <= M) { 2768865f1eaSKarl Rupp for (i = 0; i < sDist; i++) idx[nn++] = xe + i; 2778865f1eaSKarl Rupp } else { 27845b6f7e9SBarry Smith for (i = xe; i < M; i++) idx[nn++] = i; 2798865f1eaSKarl Rupp } 28047c6ae99SBarry Smith } 28147c6ae99SBarry Smith 2829566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, nn - IXs + Xs, &idx[IXs - Xs], PETSC_USE_POINTER, &from)); 2839566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global, from, local, to, >ol)); 2849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 2859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 2869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 2879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 28847c6ae99SBarry Smith 2899371c9d4SSatish Balay dd->xs = dof * xs; 2909371c9d4SSatish Balay dd->xe = dof * xe; 2919371c9d4SSatish Balay dd->ys = 0; 2929371c9d4SSatish Balay dd->ye = 1; 2939371c9d4SSatish Balay dd->zs = 0; 2949371c9d4SSatish Balay dd->ze = 1; 2959371c9d4SSatish Balay dd->Xs = dof * Xs; 2969371c9d4SSatish Balay dd->Xe = dof * Xe; 2979371c9d4SSatish Balay dd->Ys = 0; 2989371c9d4SSatish Balay dd->Ye = 1; 2999371c9d4SSatish Balay dd->Zs = 0; 3009371c9d4SSatish Balay dd->Ze = 1; 30147c6ae99SBarry Smith 30247c6ae99SBarry Smith dd->gtol = gtol; 303d6e23781SBarry Smith dd->base = dof * xs; 3049a42bb27SBarry Smith da->ops->view = DMView_DA_1d; 30547c6ae99SBarry Smith 30647c6ae99SBarry Smith /* 30747c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 30847c6ae99SBarry Smith of VecSetValuesLocal(). 30947c6ae99SBarry Smith */ 3108865f1eaSKarl Rupp for (i = 0; i < Xe - IXe; i++) idx[nn++] = -1; /* pad with -1s if needed for ghosted case*/ 311ce00eea3SSatish Balay 3129566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, dof, nn, idx, PETSC_OWN_POINTER, &da->ltogmap)); 31347c6ae99SBarry Smith 31447c6ae99SBarry Smith PetscFunctionReturn(0); 31547c6ae99SBarry Smith } 316d7bf68aeSBarry Smith 31747c6ae99SBarry Smith /*@C 318aa219208SBarry Smith DMDACreate1d - Creates an object that will manage the communication of one-dimensional 31947c6ae99SBarry Smith regular array data that is distributed across some processors. 32047c6ae99SBarry Smith 321d083f849SBarry Smith Collective 32247c6ae99SBarry Smith 32347c6ae99SBarry Smith Input Parameters: 32447c6ae99SBarry Smith + comm - MPI communicator 3251321219cSEthan Coon . bx - type of ghost cells at the boundary the array should have, if any. Use 326*dce8aebaSBarry Smith `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, or `DM_BOUNDARY_PERIODIC`. 327e43dc8daSBarry Smith . M - global dimension of the array (that is the number of grid points) 32847c6ae99SBarry Smith from the command line with -da_grid_x <M>) 32947c6ae99SBarry Smith . dof - number of degrees of freedom per node 33047c6ae99SBarry Smith . s - stencil width 33147c6ae99SBarry Smith - lx - array containing number of nodes in the X direction on each processor, 3320298fd71SBarry Smith or NULL. If non-null, must be of length as the number of processes in the MPI_Comm. 333e43dc8daSBarry Smith The sum of these entries must equal M 33447c6ae99SBarry Smith 33547c6ae99SBarry Smith Output Parameter: 33647c6ae99SBarry Smith . da - the resulting distributed array object 33747c6ae99SBarry Smith 338*dce8aebaSBarry Smith Options Database Keys: 339*dce8aebaSBarry Smith + -dm_view - Calls `DMView()` at the conclusion of `DMDACreate1d()` 340e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 341e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement factor 342*dce8aebaSBarry Smith - -da_refine <n> - refine the `DMDA` n times before creating it 34347c6ae99SBarry Smith 34447c6ae99SBarry Smith Level: beginner 34547c6ae99SBarry Smith 34647c6ae99SBarry Smith Notes: 347*dce8aebaSBarry Smith The array data itself is NOT stored in the `DMDA`, it is stored in `Vec` objects; 348*dce8aebaSBarry Smith The appropriate vector objects can be obtained with calls to `DMCreateGlobalVector()` 349*dce8aebaSBarry Smith and `DMCreateLocalVector()` and calls to `VecDuplicate()` if more are needed. 35047c6ae99SBarry Smith 351*dce8aebaSBarry Smith You must call `DMSetUp()` after this call before using this `DM`. 352897f7067SBarry Smith 353*dce8aebaSBarry Smith If you wish to use the options database to change values in the `DMDA` call `DMSetFromOptions()` after this call 354*dce8aebaSBarry Smith but before `DMSetUp()`. 355897f7067SBarry Smith 356*dce8aebaSBarry Smith .seealso: `DMDA`, `DM`, `DMDestroy()`, `DMView()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMGlobalToLocalBegin()`, `DMDASetRefinementFactor()`, 357db781477SPatrick Sanan `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDAGetRefinementFactor()`, 358db781477SPatrick Sanan `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`, 359db781477SPatrick Sanan `DMStagCreate1d()` 36047c6ae99SBarry Smith @*/ 361d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreate1d(MPI_Comm comm, DMBoundaryType bx, PetscInt M, PetscInt dof, PetscInt s, const PetscInt lx[], DM *da) 362d71ae5a4SJacob Faibussowitsch { 36347c6ae99SBarry Smith PetscMPIInt size; 36447c6ae99SBarry Smith 36547c6ae99SBarry Smith PetscFunctionBegin; 3669566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 3679566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 1)); 3689566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, 1, 1)); 3699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 3709566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, size, PETSC_DECIDE, PETSC_DECIDE)); 3719566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE)); 3729566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 3739566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 3749566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, NULL, NULL)); 37547c6ae99SBarry Smith PetscFunctionReturn(0); 37647c6ae99SBarry Smith } 377