147c6ae99SBarry Smith /* 247c6ae99SBarry Smith Code for manipulating distributed regular 1d arrays in parallel. 347c6ae99SBarry Smith This file was created by Peter Mell 6/30/95 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 747c6ae99SBarry Smith 89804daf3SBarry Smith #include <petscdraw.h> 9d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_DA_1d(DM da, PetscViewer viewer) 10d71ae5a4SJacob Faibussowitsch { 1147c6ae99SBarry Smith PetscMPIInt rank; 12*9f196a02SMartin Diehl PetscBool isascii, isdraw, isglvis, isbinary; 1347c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 14d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 159a42bb27SBarry Smith PetscBool ismatlab; 169a42bb27SBarry Smith #endif 1747c6ae99SBarry Smith 1847c6ae99SBarry Smith PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da), &rank)); 2047c6ae99SBarry Smith 21*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis)); 249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 25d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERMATLAB, &ismatlab)); 279a42bb27SBarry Smith #endif 28*9f196a02SMartin Diehl if (isascii) { 2947c6ae99SBarry Smith PetscViewerFormat format; 3047c6ae99SBarry Smith 319566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 3276a8abe0SBarry Smith if (format == PETSC_VIEWER_LOAD_BALANCE) { 331690c2aeSBarry Smith PetscInt i, nmax = 0, nmin = PETSC_INT_MAX, navg = 0, *nz, nzlocal; 3476a8abe0SBarry Smith DMDALocalInfo info; 3576a8abe0SBarry Smith PetscMPIInt size; 369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da), &size)); 379566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 3876a8abe0SBarry Smith nzlocal = info.xm; 399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &nz)); 409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nzlocal, 1, MPIU_INT, nz, 1, MPIU_INT, PetscObjectComm((PetscObject)da))); 41835f2295SStefano Zampini for (i = 0; i < size; i++) { 4276a8abe0SBarry Smith nmax = PetscMax(nmax, nz[i]); 4376a8abe0SBarry Smith nmin = PetscMin(nmin, nz[i]); 4476a8abe0SBarry Smith navg += nz[i]; 4576a8abe0SBarry Smith } 469566063dSJacob Faibussowitsch PetscCall(PetscFree(nz)); 4776a8abe0SBarry Smith navg = navg / size; 4863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Load Balance - Grid Points: Min %" PetscInt_FMT " avg %" PetscInt_FMT " max %" PetscInt_FMT "\n", nmin, navg, nmax)); 493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5076a8abe0SBarry Smith } 51ce78bad3SBarry Smith if (format != PETSC_VIEWER_ASCII_GLVIS) { 52aa219208SBarry Smith DMDALocalInfo info; 539566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 5563a3b9bcSJacob 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)); 5663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "X range of indices: %" PetscInt_FMT " %" PetscInt_FMT "\n", info.xs, info.xs + info.xm)); 579566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 591baa6e33SBarry Smith } else if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMView_DA_GLVis(da, viewer)); 6047c6ae99SBarry Smith } else if (isdraw) { 6147c6ae99SBarry Smith PetscDraw draw; 6247c6ae99SBarry Smith double ymin = -1, ymax = 1, xmin = -1, xmax = dd->M, x; 6347c6ae99SBarry Smith PetscInt base; 6447c6ae99SBarry Smith char node[10]; 6547c6ae99SBarry Smith PetscBool isnull; 6647c6ae99SBarry Smith 679566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 689566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw, &isnull)); 693ba16761SJacob Faibussowitsch if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 705b399a63SLisandro Dalcin 719566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 729566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 739566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw, xmin, ymin, xmax, ymax)); 7447c6ae99SBarry Smith 75d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 7647c6ae99SBarry Smith /* first processor draws all node lines */ 77dd400576SPatrick Sanan if (rank == 0) { 7847c6ae99SBarry Smith PetscInt xmin_tmp; 799371c9d4SSatish Balay ymin = 0.0; 809371c9d4SSatish Balay ymax = 0.3; 8148a46eb9SPierre 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)); 829371c9d4SSatish Balay xmin = 0.0; 839371c9d4SSatish Balay xmax = dd->M - 1; 849566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ymin, PETSC_DRAW_BLACK)); 859566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymax, xmax, ymax, PETSC_DRAW_BLACK)); 8647c6ae99SBarry Smith } 87d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 889566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 899566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 9047c6ae99SBarry Smith 91d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 9247c6ae99SBarry Smith /* draw my box */ 939371c9d4SSatish Balay ymin = 0; 949371c9d4SSatish Balay ymax = 0.3; 959371c9d4SSatish Balay xmin = dd->xs / dd->w; 969371c9d4SSatish Balay xmax = (dd->xe / dd->w) - 1; 979566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ymin, PETSC_DRAW_RED)); 989566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymin, xmin, ymax, PETSC_DRAW_RED)); 999566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymax, xmax, ymax, PETSC_DRAW_RED)); 1009566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmax, ymin, xmax, ymax, PETSC_DRAW_RED)); 10147c6ae99SBarry Smith /* Put in index numbers */ 10247c6ae99SBarry Smith base = dd->base / dd->w; 10347c6ae99SBarry Smith for (x = xmin; x <= xmax; x++) { 104835f2295SStefano Zampini PetscCall(PetscSNPrintf(node, sizeof(node), "%" PetscInt_FMT, base++)); 1059566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, x, ymin, PETSC_DRAW_RED, node)); 10647c6ae99SBarry Smith } 107d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1089566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1099566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1109566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 111c9493c35SStefano Zampini } else if (isglvis) { 1129566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da, viewer)); 1139a42bb27SBarry Smith } else if (isbinary) { 1149566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da, viewer)); 115d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 1169a42bb27SBarry Smith } else if (ismatlab) { 1179566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da, viewer)); 1189a42bb27SBarry Smith #endif 11911aeaf0aSBarry Smith } 1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12147c6ae99SBarry Smith } 12247c6ae99SBarry Smith 123d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUp_DA_1D(DM da) 124d71ae5a4SJacob Faibussowitsch { 12547c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 12647c6ae99SBarry Smith const PetscInt M = dd->M; 12747c6ae99SBarry Smith const PetscInt dof = dd->w; 12847c6ae99SBarry Smith const PetscInt s = dd->s; 12945b6f7e9SBarry Smith const PetscInt sDist = s; /* stencil distance in points */ 13047c6ae99SBarry Smith const PetscInt *lx = dd->lx; 131bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 13247c6ae99SBarry Smith MPI_Comm comm; 13347c6ae99SBarry Smith Vec local, global; 134bd1fc5aeSBarry Smith VecScatter gtol; 13547c6ae99SBarry Smith IS to, from; 13647c6ae99SBarry Smith PetscBool flg1 = PETSC_FALSE, flg2 = PETSC_FALSE; 13747c6ae99SBarry Smith PetscMPIInt rank, size; 138bd1fc5aeSBarry Smith PetscInt i, *idx, nn, left, xs, xe, x, Xs, Xe, start, m, IXs, IXe; 13947c6ae99SBarry Smith 14047c6ae99SBarry Smith PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 1429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 14447c6ae99SBarry Smith 1457d310018SBarry Smith dd->p = 1; 1467d310018SBarry Smith dd->n = 1; 14747c6ae99SBarry Smith dd->m = size; 14847c6ae99SBarry Smith m = dd->m; 14947c6ae99SBarry Smith 15047c6ae99SBarry Smith if (s > 0) { 15147c6ae99SBarry Smith /* if not communicating data then should be ok to have nothing on some processes */ 15263a3b9bcSJacob Faibussowitsch PetscCheck(M >= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "More processes than data points! %" PetscInt_FMT " %" PetscInt_FMT, m, M); 15363a3b9bcSJacob 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); 15447c6ae99SBarry Smith } 15547c6ae99SBarry Smith 15647c6ae99SBarry Smith /* 15747c6ae99SBarry Smith Determine locally owned region 15847c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 15947c6ae99SBarry Smith */ 16047c6ae99SBarry Smith if (!lx) { 1619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 1629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_partition_blockcomm", &flg1, NULL)); 1639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_partition_nodes_at_end", &flg2, NULL)); 16447c6ae99SBarry Smith if (flg1) { /* Block Comm type Distribution */ 16547c6ae99SBarry Smith xs = rank * M / m; 16647c6ae99SBarry Smith x = (rank + 1) * M / m - xs; 16747c6ae99SBarry Smith } else if (flg2) { /* The odd nodes are evenly distributed across last nodes */ 16847c6ae99SBarry Smith x = (M + rank) / m; 1698865f1eaSKarl Rupp if (M / m == x) xs = rank * x; 1708865f1eaSKarl Rupp else xs = rank * (x - 1) + (M + rank) % (x * m); 17147c6ae99SBarry Smith } else { /* The odd nodes are evenly distributed across the first k nodes */ 17247c6ae99SBarry Smith /* Regular PETSc Distribution */ 17347c6ae99SBarry Smith x = M / m + ((M % m) > rank); 174835f2295SStefano Zampini if (rank >= (M % m)) xs = (rank * (M / m) + M % m); 175835f2295SStefano Zampini else xs = rank * (M / m) + rank; 17647c6ae99SBarry Smith } 1779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&xs, 1, MPIU_INT, dd->lx, 1, MPIU_INT, comm)); 178fe16a2e9SBarry Smith for (i = 0; i < m - 1; i++) dd->lx[i] = dd->lx[i + 1] - dd->lx[i]; 179fe16a2e9SBarry Smith dd->lx[m - 1] = M - dd->lx[m - 1]; 18047c6ae99SBarry Smith } else { 18147c6ae99SBarry Smith x = lx[rank]; 18247c6ae99SBarry Smith xs = 0; 1838865f1eaSKarl Rupp for (i = 0; i < rank; i++) xs += lx[i]; 18447c6ae99SBarry Smith /* verify that data user provided is consistent */ 18547c6ae99SBarry Smith left = xs; 1868865f1eaSKarl Rupp for (i = rank; i < size; i++) left += lx[i]; 18763a3b9bcSJacob 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); 18847c6ae99SBarry Smith } 18947c6ae99SBarry Smith 190bcea557cSEthan Coon /* 191bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 192bcea557cSEthan Coon the domain more than once 193bcea557cSEthan Coon */ 19463a3b9bcSJacob 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); 195bcea557cSEthan Coon 19647c6ae99SBarry Smith xe = xs + x; 19747c6ae99SBarry Smith 19888661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 199d9c9ebe5SPeter Brune if (xs - sDist > 0) { 200d9c9ebe5SPeter Brune Xs = xs - sDist; 201d9c9ebe5SPeter Brune IXs = xs - sDist; 20288661749SPeter Brune } else { 2038865f1eaSKarl Rupp if (bx) Xs = xs - sDist; 2048865f1eaSKarl Rupp else Xs = 0; 20588661749SPeter Brune IXs = 0; 20688661749SPeter Brune } 20745b6f7e9SBarry Smith if (xe + sDist <= M) { 208d9c9ebe5SPeter Brune Xe = xe + sDist; 209d9c9ebe5SPeter Brune IXe = xe + sDist; 21088661749SPeter Brune } else { 2118865f1eaSKarl Rupp if (bx) Xe = xe + sDist; 21245b6f7e9SBarry Smith else Xe = M; 21345b6f7e9SBarry Smith IXe = M; 21447c6ae99SBarry Smith } 21547c6ae99SBarry Smith 216bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 217d9c9ebe5SPeter Brune Xs = xs - sDist; 218d9c9ebe5SPeter Brune Xe = xe + sDist; 219d9c9ebe5SPeter Brune IXs = xs - sDist; 220d9c9ebe5SPeter Brune IXe = xe + sDist; 221ce00eea3SSatish Balay } 222ce00eea3SSatish Balay 22347c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 22445b6f7e9SBarry Smith dd->Nlocal = dof * x; 2259566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, dof, dd->Nlocal, PETSC_DECIDE, NULL, &global)); 22645b6f7e9SBarry Smith dd->nlocal = dof * (Xe - Xs); 2279566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dof, dd->nlocal, NULL, &local)); 22847c6ae99SBarry Smith 2299566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(global, &start, NULL)); 23047c6ae99SBarry Smith 23147c6ae99SBarry Smith /* Create Global to Local Vector Scatter Context */ 23247c6ae99SBarry Smith /* global to local must retrieve ghost points */ 2339566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, dof * (IXe - IXs), dof * (IXs - Xs), 1, &to)); 23447c6ae99SBarry Smith 2359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(x + 2 * sDist, &idx)); 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, >ol)); 2829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 2839566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 2849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 2859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 28647c6ae99SBarry Smith 2879371c9d4SSatish Balay dd->xs = dof * xs; 2889371c9d4SSatish Balay dd->xe = dof * xe; 2899371c9d4SSatish Balay dd->ys = 0; 2909371c9d4SSatish Balay dd->ye = 1; 2919371c9d4SSatish Balay dd->zs = 0; 2929371c9d4SSatish Balay dd->ze = 1; 2939371c9d4SSatish Balay dd->Xs = dof * Xs; 2949371c9d4SSatish Balay dd->Xe = dof * Xe; 2959371c9d4SSatish Balay dd->Ys = 0; 2969371c9d4SSatish Balay dd->Ye = 1; 2979371c9d4SSatish Balay dd->Zs = 0; 2989371c9d4SSatish Balay dd->Ze = 1; 29947c6ae99SBarry Smith 30047c6ae99SBarry Smith dd->gtol = gtol; 301d6e23781SBarry Smith dd->base = dof * xs; 3029a42bb27SBarry Smith da->ops->view = DMView_DA_1d; 30347c6ae99SBarry Smith 30447c6ae99SBarry Smith /* 30547c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 30647c6ae99SBarry Smith of VecSetValuesLocal(). 30747c6ae99SBarry Smith */ 3088865f1eaSKarl Rupp for (i = 0; i < Xe - IXe; i++) idx[nn++] = -1; /* pad with -1s if needed for ghosted case*/ 309ce00eea3SSatish Balay 3109566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, dof, nn, idx, PETSC_OWN_POINTER, &da->ltogmap)); 3113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31247c6ae99SBarry Smith } 313d7bf68aeSBarry Smith 314cc4c1da9SBarry Smith /*@ 315aa219208SBarry Smith DMDACreate1d - Creates an object that will manage the communication of one-dimensional 31612b4a537SBarry Smith regular array data that is distributed across one or mpre MPI processes. 31747c6ae99SBarry Smith 318d083f849SBarry 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 323dce8aebaSBarry Smith `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 . dof - number of degrees of freedom per node 32647c6ae99SBarry Smith . s - stencil width 32747c6ae99SBarry Smith - lx - array containing number of nodes in the X direction on each processor, 32812b4a537SBarry Smith or `NULL`. If non-null, must be of length as the number of processes in the MPI_Comm. 32912b4a537SBarry Smith The sum of these entries must equal `M` 33047c6ae99SBarry Smith 33147c6ae99SBarry Smith Output Parameter: 33247c6ae99SBarry Smith . da - the resulting distributed array object 33347c6ae99SBarry Smith 334dce8aebaSBarry Smith Options Database Keys: 335dce8aebaSBarry Smith + -dm_view - Calls `DMView()` at the conclusion of `DMDACreate1d()` 336e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 337e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement factor 338dce8aebaSBarry Smith - -da_refine <n> - refine the `DMDA` n times before creating it 33947c6ae99SBarry Smith 34047c6ae99SBarry Smith Level: beginner 34147c6ae99SBarry Smith 34247c6ae99SBarry Smith Notes: 343dce8aebaSBarry Smith The array data itself is NOT stored in the `DMDA`, it is stored in `Vec` objects; 344dce8aebaSBarry Smith The appropriate vector objects can be obtained with calls to `DMCreateGlobalVector()` 345dce8aebaSBarry Smith and `DMCreateLocalVector()` and calls to `VecDuplicate()` if more are needed. 34647c6ae99SBarry Smith 347dce8aebaSBarry Smith You must call `DMSetUp()` after this call before using this `DM`. 348897f7067SBarry Smith 349dce8aebaSBarry Smith If you wish to use the options database to change values in the `DMDA` call `DMSetFromOptions()` after this call 350dce8aebaSBarry Smith but before `DMSetUp()`. 351897f7067SBarry Smith 35212b4a537SBarry Smith .seealso: [](sec_struct), `DMDA`, `DM`, `DMDestroy()`, `DMView()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMGlobalToLocalBegin()`, `DMDASetRefinementFactor()`, 353db781477SPatrick Sanan `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDAGetRefinementFactor()`, 354db781477SPatrick Sanan `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`, 35512b4a537SBarry Smith `DMStagCreate1d()`, `DMBoundaryType` 35647c6ae99SBarry Smith @*/ 357d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreate1d(MPI_Comm comm, DMBoundaryType bx, PetscInt M, PetscInt dof, PetscInt s, const PetscInt lx[], DM *da) 358d71ae5a4SJacob Faibussowitsch { 35947c6ae99SBarry Smith PetscMPIInt size; 36047c6ae99SBarry Smith 36147c6ae99SBarry Smith PetscFunctionBegin; 3629566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 3639566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 1)); 3649566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, 1, 1)); 3659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 3669566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, size, PETSC_DECIDE, PETSC_DECIDE)); 3679566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE)); 3689566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 3699566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 3709566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, NULL, NULL)); 3713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37247c6ae99SBarry Smith } 373