Actual source code: daview.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  2: /*
  3:   Code for manipulating distributed regular arrays in parallel.
  4: */

  6: #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/

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

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

 24:   MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);
 25:   if (!rank) {
 26:     DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&dof,&swidth,&bx,&by,&bz,&stencil);
 27:     mx   = mxCreateStructMatrix(1,1,8,(const char**)fnames);
 28:     if (!mx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unable to generate MATLAB struct array to hold DMDA informations");
 29:     mxSetFieldByNumber(mx,0,0,mxCreateDoubleScalar((double)dim));
 30:     mxSetFieldByNumber(mx,0,1,mxCreateDoubleScalar((double)m));
 31:     mxSetFieldByNumber(mx,0,2,mxCreateDoubleScalar((double)n));
 32:     mxSetFieldByNumber(mx,0,3,mxCreateDoubleScalar((double)p));
 33:     mxSetFieldByNumber(mx,0,4,mxCreateDoubleScalar((double)dof));
 34:     mxSetFieldByNumber(mx,0,5,mxCreateDoubleScalar((double)swidth));
 35:     mxSetFieldByNumber(mx,0,6,mxCreateDoubleScalar((double)bx));
 36:     mxSetFieldByNumber(mx,0,7,mxCreateDoubleScalar((double)by));
 37:     mxSetFieldByNumber(mx,0,8,mxCreateDoubleScalar((double)bz));
 38:     mxSetFieldByNumber(mx,0,9,mxCreateDoubleScalar((double)stencil));
 39:     PetscObjectName((PetscObject)da);
 40:     PetscViewerMatlabPutVariable(viewer,((PetscObject)da)->name,mx);
 41:   }
 42:   return(0);
 43: }
 44: #endif

 48: PetscErrorCode DMView_DA_Binary(DM da,PetscViewer viewer)
 49: {
 50:   PetscErrorCode   ierr;
 51:   PetscMPIInt      rank;
 52:   PetscInt         dim,m,n,p,dof,swidth,M,N,P;
 53:   DMDAStencilType  stencil;
 54:   DMBoundaryType   bx,by,bz;
 55:   MPI_Comm         comm;
 56:   PetscBool        coors = PETSC_FALSE;

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

 61:   DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&dof,&swidth,&bx,&by,&bz,&stencil);
 62:   MPI_Comm_rank(comm,&rank);
 63:   if (!rank) {
 64:     PetscViewerBinaryWrite(viewer,&dim,1,PETSC_INT,PETSC_FALSE);
 65:     PetscViewerBinaryWrite(viewer,&m,1,PETSC_INT,PETSC_FALSE);
 66:     PetscViewerBinaryWrite(viewer,&n,1,PETSC_INT,PETSC_FALSE);
 67:     PetscViewerBinaryWrite(viewer,&p,1,PETSC_INT,PETSC_FALSE);
 68:     PetscViewerBinaryWrite(viewer,&dof,1,PETSC_INT,PETSC_FALSE);
 69:     PetscViewerBinaryWrite(viewer,&swidth,1,PETSC_INT,PETSC_FALSE);
 70:     PetscViewerBinaryWrite(viewer,&bx,1,PETSC_ENUM,PETSC_FALSE);
 71:     PetscViewerBinaryWrite(viewer,&by,1,PETSC_ENUM,PETSC_FALSE);
 72:     PetscViewerBinaryWrite(viewer,&bz,1,PETSC_ENUM,PETSC_FALSE);
 73:     PetscViewerBinaryWrite(viewer,&stencil,1,PETSC_ENUM,PETSC_FALSE);
 74:     if (da->coordinates) coors = PETSC_TRUE;
 75:     PetscViewerBinaryWrite(viewer,&coors,1,PETSC_BOOL,PETSC_FALSE);
 76:   }

 78:   /* save the coordinates if they exist to disk (in the natural ordering) */
 79:   if (da->coordinates) {
 80:     VecView(da->coordinates,viewer);
 81:   }
 82:   return(0);
 83: }

 87: PetscErrorCode DMView_DA_VTK(DM da, PetscViewer viewer)
 88: {
 89:   PetscInt       dim, dof, M = 0, N = 0, P = 0;

 93:   DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL);
 94:   if (!da->coordinates) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP, "VTK output requires DMDA coordinates.");
 95:   /* Write Header */
 96:   PetscViewerASCIIPrintf(viewer,"# vtk DataFile Version 2.0\n");
 97:   PetscViewerASCIIPrintf(viewer,"Structured Mesh Example\n");
 98:   PetscViewerASCIIPrintf(viewer,"ASCII\n");
 99:   PetscViewerASCIIPrintf(viewer,"DATASET STRUCTURED_GRID\n");
100:   PetscViewerASCIIPrintf(viewer,"DIMENSIONS %d %d %d\n", M, N, P);
101:   PetscViewerASCIIPrintf(viewer,"POINTS %d double\n", M*N*P);
102:   if (da->coordinates) {
103:     DM  dac;
104:     Vec natural;

106:     DMGetCoordinateDM(da, &dac);
107:     DMDACreateNaturalVector(dac, &natural);
108:     PetscObjectSetOptionsPrefix((PetscObject) natural, "coor_");
109:     DMDAGlobalToNaturalBegin(dac, da->coordinates, INSERT_VALUES, natural);
110:     DMDAGlobalToNaturalEnd(dac, da->coordinates, INSERT_VALUES, natural);
111:     PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_VTK_COORDS);
112:     VecView(natural, viewer);
113:     PetscViewerPopFormat(viewer);
114:     VecDestroy(&natural);
115:   }
116:   return(0);
117: }

121: /*@C
122:    DMDAGetInfo - Gets information about a given distributed array.

124:    Not Collective

126:    Input Parameter:
127: .  da - the distributed array

129:    Output Parameters:
130: +  dim      - dimension of the distributed array (1, 2, or 3)
131: .  M, N, P  - global dimension in each direction of the array
132: .  m, n, p  - corresponding number of procs in each dimension
133: .  dof      - number of degrees of freedom per node
134: .  s        - stencil width
135: .  bx,by,bz - type of ghost nodes at boundary, one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED,
136:               DM_BOUNDARY_MIRROR, DM_BOUNDARY_PERIODIC
137: -  st       - stencil type, either DMDA_STENCIL_STAR or DMDA_STENCIL_BOX

139:    Level: beginner

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

144: .keywords: distributed array, get, information

146: .seealso: DMView(), DMDAGetCorners(), DMDAGetLocalInfo()
147: @*/
148: 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)
149: {
150:   DM_DA *dd = (DM_DA*)da->data;

154:   if (dim) *dim = da->dim;
155:   if (M) {
156:     if (dd->Mo < 0) *M = dd->M;
157:     else *M = dd->Mo;
158:   }
159:   if (N) {
160:     if (dd->No < 0) *N = dd->N;
161:     else *N = dd->No;
162:   }
163:   if (P) {
164:     if (dd->Po < 0) *P = dd->P;
165:     else *P = dd->Po;
166:   }
167:   if (m) *m = dd->m;
168:   if (n) *n = dd->n;
169:   if (p) *p = dd->p;
170:   if (dof) *dof = dd->w;
171:   if (s) *s = dd->s;
172:   if (bx) *bx = dd->bx;
173:   if (by) *by = dd->by;
174:   if (bz) *bz = dd->bz;
175:   if (st) *st = dd->stencil_type;
176:   return(0);
177: }

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

184:    Not Collective

186:    Input Parameter:
187: .  da - the distributed array

189:    Output Parameters:
190: .  dainfo - structure containing the information

192:    Level: beginner

194: .keywords: distributed array, get, information

196: .seealso: DMDAGetInfo(), DMDAGetCorners()
197: @*/
198: PetscErrorCode  DMDAGetLocalInfo(DM da,DMDALocalInfo *info)
199: {
200:   PetscInt w;
201:   DM_DA    *dd = (DM_DA*)da->data;

206:   info->da  = da;
207:   info->dim = da->dim;
208:   if (dd->Mo < 0) info->mx = dd->M;
209:   else info->mx = dd->Mo;
210:   if (dd->No < 0) info->my = dd->N;
211:   else info->my = dd->No;
212:   if (dd->Po < 0) info->mz = dd->P;
213:   else info->mz = dd->Po;
214:   info->dof = dd->w;
215:   info->sw  = dd->s;
216:   info->bx  = dd->bx;
217:   info->by  = dd->by;
218:   info->bz  = dd->bz;
219:   info->st  = dd->stencil_type;

221:   /* since the xs, xe ... have all been multiplied by the number of degrees
222:      of freedom per cell, w = dd->w, we divide that out before returning.*/
223:   w        = dd->w;
224:   info->xs = dd->xs/w + dd->xo;
225:   info->xm = (dd->xe - dd->xs)/w;
226:   /* the y and z have NOT been multiplied by w */
227:   info->ys = dd->ys + dd->yo;
228:   info->ym = (dd->ye - dd->ys);
229:   info->zs = dd->zs + dd->zo;
230:   info->zm = (dd->ze - dd->zs);

232:   info->gxs = dd->Xs/w + dd->xo;
233:   info->gxm = (dd->Xe - dd->Xs)/w;
234:   /* the y and z have NOT been multiplied by w */
235:   info->gys = dd->Ys + dd->yo;
236:   info->gym = (dd->Ye - dd->Ys);
237:   info->gzs = dd->Zs + dd->zo;
238:   info->gzm = (dd->Ze - dd->Zs);
239:   return(0);
240: }