Actual source code: daview.c


  2: /*
  3:   Code for manipulating distributed regular arrays in parallel.
  4: */

  6: #include <petsc/private/dmdaimpl.h>

  8: #if defined(PETSC_HAVE_MATLAB)
  9:   #include <mat.h> /* MATLAB include file */

 11: PetscErrorCode DMView_DA_Matlab(DM da, PetscViewer viewer)
 12: {
 13:   PetscMPIInt     rank;
 14:   PetscInt        dim, m, n, p, dof, swidth;
 15:   DMDAStencilType stencil;
 16:   DMBoundaryType  bx, by, bz;
 17:   mxArray        *mx;
 18:   const char     *fnames[] = {"dimension", "m", "n", "p", "dof", "stencil_width", "bx", "by", "bz", "stencil_type"};

 20:   MPI_Comm_rank(PetscObjectComm((PetscObject)da), &rank);
 21:   if (rank == 0) {
 22:     DMDAGetInfo(da, &dim, &m, &n, &p, 0, 0, 0, &dof, &swidth, &bx, &by, &bz, &stencil);
 23:     mx = mxCreateStructMatrix(1, 1, 8, (const char **)fnames);
 25:     mxSetFieldByNumber(mx, 0, 0, mxCreateDoubleScalar((double)dim));
 26:     mxSetFieldByNumber(mx, 0, 1, mxCreateDoubleScalar((double)m));
 27:     mxSetFieldByNumber(mx, 0, 2, mxCreateDoubleScalar((double)n));
 28:     mxSetFieldByNumber(mx, 0, 3, mxCreateDoubleScalar((double)p));
 29:     mxSetFieldByNumber(mx, 0, 4, mxCreateDoubleScalar((double)dof));
 30:     mxSetFieldByNumber(mx, 0, 5, mxCreateDoubleScalar((double)swidth));
 31:     mxSetFieldByNumber(mx, 0, 6, mxCreateDoubleScalar((double)bx));
 32:     mxSetFieldByNumber(mx, 0, 7, mxCreateDoubleScalar((double)by));
 33:     mxSetFieldByNumber(mx, 0, 8, mxCreateDoubleScalar((double)bz));
 34:     mxSetFieldByNumber(mx, 0, 9, mxCreateDoubleScalar((double)stencil));
 35:     PetscObjectName((PetscObject)da);
 36:     PetscViewerMatlabPutVariable(viewer, ((PetscObject)da)->name, mx);
 37:   }
 38:   return 0;
 39: }
 40: #endif

 42: PetscErrorCode DMView_DA_Binary(DM da, PetscViewer viewer)
 43: {
 44:   PetscMPIInt     rank;
 45:   PetscInt        dim, m, n, p, dof, swidth, M, N, P;
 46:   DMDAStencilType stencil;
 47:   DMBoundaryType  bx, by, bz;
 48:   MPI_Comm        comm;
 49:   PetscBool       coors = PETSC_FALSE;
 50:   Vec             coordinates;

 52:   PetscObjectGetComm((PetscObject)da, &comm);

 54:   DMDAGetInfo(da, &dim, &m, &n, &p, &M, &N, &P, &dof, &swidth, &bx, &by, &bz, &stencil);
 55:   MPI_Comm_rank(comm, &rank);
 56:   DMGetCoordinates(da, &coordinates);
 57:   if (rank == 0) {
 58:     PetscViewerBinaryWrite(viewer, &dim, 1, PETSC_INT);
 59:     PetscViewerBinaryWrite(viewer, &m, 1, PETSC_INT);
 60:     PetscViewerBinaryWrite(viewer, &n, 1, PETSC_INT);
 61:     PetscViewerBinaryWrite(viewer, &p, 1, PETSC_INT);
 62:     PetscViewerBinaryWrite(viewer, &dof, 1, PETSC_INT);
 63:     PetscViewerBinaryWrite(viewer, &swidth, 1, PETSC_INT);
 64:     PetscViewerBinaryWrite(viewer, &bx, 1, PETSC_ENUM);
 65:     PetscViewerBinaryWrite(viewer, &by, 1, PETSC_ENUM);
 66:     PetscViewerBinaryWrite(viewer, &bz, 1, PETSC_ENUM);
 67:     PetscViewerBinaryWrite(viewer, &stencil, 1, PETSC_ENUM);
 68:     if (coordinates) coors = PETSC_TRUE;
 69:     PetscViewerBinaryWrite(viewer, &coors, 1, PETSC_BOOL);
 70:   }

 72:   /* save the coordinates if they exist to disk (in the natural ordering) */
 73:   if (coordinates) VecView(coordinates, viewer);
 74:   return 0;
 75: }

 77: PetscErrorCode DMView_DA_VTK(DM da, PetscViewer viewer)
 78: {
 79:   Vec      coordinates;
 80:   PetscInt dim, dof, M = 0, N = 0, P = 0;

 82:   DMGetCoordinates(da, &coordinates);
 83:   DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL);
 85:   /* Write Header */
 86:   PetscViewerASCIIPrintf(viewer, "# vtk DataFile Version 2.0\n");
 87:   PetscViewerASCIIPrintf(viewer, "Structured Mesh Example\n");
 88:   PetscViewerASCIIPrintf(viewer, "ASCII\n");
 89:   PetscViewerASCIIPrintf(viewer, "DATASET STRUCTURED_GRID\n");
 90:   PetscViewerASCIIPrintf(viewer, "DIMENSIONS %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", M, N, P);
 91:   PetscViewerASCIIPrintf(viewer, "POINTS %" PetscInt_FMT " double\n", M * N * P);
 92:   if (coordinates) {
 93:     DM  dac;
 94:     Vec natural;

 96:     DMGetCoordinateDM(da, &dac);
 97:     DMDACreateNaturalVector(dac, &natural);
 98:     PetscObjectSetOptionsPrefix((PetscObject)natural, "coor_");
 99:     DMDAGlobalToNaturalBegin(dac, coordinates, INSERT_VALUES, natural);
100:     DMDAGlobalToNaturalEnd(dac, coordinates, INSERT_VALUES, natural);
101:     PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_VTK_COORDS_DEPRECATED);
102:     VecView(natural, viewer);
103:     PetscViewerPopFormat(viewer);
104:     VecDestroy(&natural);
105:   }
106:   return 0;
107: }

109: /*@C
110:    DMDAGetInfo - Gets information about a given distributed array.

112:    Not Collective

114:    Input Parameter:
115: .  da - the distributed array

117:    Output Parameters:
118: +  dim      - dimension of the distributed array (1, 2, or 3)
119: .  M        - global dimension in first direction of the array
120: .  N        - global dimension in second direction of the array
121: .  P        - global dimension in third direction of the array
122: .  m        - corresponding number of procs in first dimension
123: .  n        - corresponding number of procs in second dimension
124: .  p        - corresponding number of procs in third dimension
125: .  dof      - number of degrees of freedom per node
126: .  s        - stencil width
127: .  bx       - type of ghost nodes at boundary in first dimension
128: .  by       - type of ghost nodes at boundary in second dimension
129: .  bz       - type of ghost nodes at boundary in third dimension
130: -  st       - stencil type, either `DMDA_STENCIL_STAR` or `DMDA_STENCIL_BOX`

132:    Level: beginner

134:    Note:
135:    Use NULL (NULL_INTEGER in Fortran) in place of any output parameter that is not of interest.

137: .seealso: `DM`, `DMDA`, `DMView()`, `DMDAGetCorners()`, `DMDAGetLocalInfo()`
138: @*/
139: PetscErrorCode DMDAGetInfo(DM da, PetscInt *dim, PetscInt *M, PetscInt *N, PetscInt *P, PetscInt *m, PetscInt *n, PetscInt *p, PetscInt *dof, PetscInt *s, DMBoundaryType *bx, DMBoundaryType *by, DMBoundaryType *bz, DMDAStencilType *st)
140: {
141:   DM_DA *dd = (DM_DA *)da->data;

144:   if (dim) *dim = da->dim;
145:   if (M) {
146:     if (dd->Mo < 0) *M = dd->M;
147:     else *M = dd->Mo;
148:   }
149:   if (N) {
150:     if (dd->No < 0) *N = dd->N;
151:     else *N = dd->No;
152:   }
153:   if (P) {
154:     if (dd->Po < 0) *P = dd->P;
155:     else *P = dd->Po;
156:   }
157:   if (m) *m = dd->m;
158:   if (n) *n = dd->n;
159:   if (p) *p = dd->p;
160:   if (dof) *dof = dd->w;
161:   if (s) *s = dd->s;
162:   if (bx) *bx = dd->bx;
163:   if (by) *by = dd->by;
164:   if (bz) *bz = dd->bz;
165:   if (st) *st = dd->stencil_type;
166:   return 0;
167: }

169: /*@C
170:    DMDAGetLocalInfo - Gets information about a given distributed array and this processors location in it

172:    Not Collective

174:    Input Parameter:
175: .  da - the distributed array

177:    Output Parameters:
178: .  dainfo - structure containing the information

180:    Level: beginner

182:    Note:
183:     See `DMDALocalInfo` for the information that is returned

185: .seealso: `DM`, `DMDA`, `DMDAGetInfo()`, `DMDAGetCorners()`, `DMDALocalInfo`
186: @*/
187: PetscErrorCode DMDAGetLocalInfo(DM da, DMDALocalInfo *info)
188: {
189:   PetscInt w;
190:   DM_DA   *dd = (DM_DA *)da->data;

194:   info->da  = da;
195:   info->dim = da->dim;
196:   if (dd->Mo < 0) info->mx = dd->M;
197:   else info->mx = dd->Mo;
198:   if (dd->No < 0) info->my = dd->N;
199:   else info->my = dd->No;
200:   if (dd->Po < 0) info->mz = dd->P;
201:   else info->mz = dd->Po;
202:   info->dof = dd->w;
203:   info->sw  = dd->s;
204:   info->bx  = dd->bx;
205:   info->by  = dd->by;
206:   info->bz  = dd->bz;
207:   info->st  = dd->stencil_type;

209:   /* since the xs, xe ... have all been multiplied by the number of degrees
210:      of freedom per cell, w = dd->w, we divide that out before returning.*/
211:   w        = dd->w;
212:   info->xs = dd->xs / w + dd->xo;
213:   info->xm = (dd->xe - dd->xs) / w;
214:   /* the y and z have NOT been multiplied by w */
215:   info->ys = dd->ys + dd->yo;
216:   info->ym = (dd->ye - dd->ys);
217:   info->zs = dd->zs + dd->zo;
218:   info->zm = (dd->ze - dd->zs);

220:   info->gxs = dd->Xs / w + dd->xo;
221:   info->gxm = (dd->Xe - dd->Xs) / w;
222:   /* the y and z have NOT been multiplied by w */
223:   info->gys = dd->Ys + dd->yo;
224:   info->gym = (dd->Ye - dd->Ys);
225:   info->gzs = dd->Zs + dd->zo;
226:   info->gzm = (dd->Ze - dd->Zs);
227:   return 0;
228: }