Actual source code: vmatlab.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <petsc-private/viewerimpl.h>
  3: #include <mat.h>

  5: /*MC
  6:    PETSCVIEWERMATLAB - A viewer that saves the variables into a MATLAB .mat file that may be read into MATLAB
  7:        with load('filename').

  9:    Level: intermediate

 11:        Note: Currently can only save PETSc vectors to .mat files, not matrices (use the PETSCVIEWERBINARY and
 12:              ${PETSC_DIR}/share/petsc/matlab/PetscBinaryRead.m to read matrices into MATLAB).

 14:              For parallel vectors obtained with DMCreateGlobalVector() or DMGetGlobalVector() the vectors are saved to
 15:              the .mat file in natural ordering. You can use DMView() to save the DMDA information to the .mat file
 16:              the fields in the MATLAB loaded da variable give the array dimensions so you can reshape the MATLAB
 17:              vector to the same multidimensional shape as it had in PETSc for plotting etc. For example,

 19: $             In your PETSc C/C++ code (assuming a two dimensional DMDA with one degree of freedom per node)
 20: $                PetscObjectSetName((PetscObject)x,"x");
 21: $                VecView(x,PETSC_VIEWER_MATLAB_WORLD);
 22: $                PetscObjectSetName((PetscObject)da,"da");
 23: $                DMView(x,PETSC_VIEWER_MATLAB_WORLD);
 24: $             Then from MATLAB
 25: $                load('matlaboutput.mat')   % matlaboutput.mat is the default filename
 26: $                xnew = zeros(da.n,da.m);
 27: $                xnew(:) = x;    % reshape one dimensional vector back to two dimensions

 29:               If you wish to put the same variable into the .mat file several times you need to give it a new
 30:               name before each call to view.

 32:               Use PetscViewerMatlabPutArray() to just put an array of doubles into the .mat file

 34: .seealso:  PETSC_VIEWER_MATLAB_(),PETSC_VIEWER_MATLAB_SELF(), PETSC_VIEWER_MATLAB_WORLD(),PetscViewerCreate(),
 35:            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERBINARY,
 36:            PETSC_ASCII_VIEWER, PetscViewerFileSetName(), PetscViewerFileSetMode()

 38: M*/

 40: typedef struct {
 41:   MATFile       *ep;
 42:   PetscMPIInt   rank;
 43:   PetscFileMode btype;
 44: } PetscViewer_Matlab;

 48: /*@C
 49:     PetscViewerMatlabPutArray - Puts an array into the MATLAB viewer.

 51:       Not collective: only processor zero saves the array

 53:     Input Parameters:
 54: +    mfile - the viewer
 55: .    m,n - the dimensions of the array
 56: .    array - the array (represented in one dimension)
 57: -    name - the name of the array

 59:    Level: advanced

 61:      Notes: Only writes array values on processor 0.

 63: @*/
 64: PetscErrorCode  PetscViewerMatlabPutArray(PetscViewer mfile,int m,int n,const PetscScalar *array,const char *name)
 65: {
 66:   PetscErrorCode     ierr;
 67:   PetscViewer_Matlab *ml = (PetscViewer_Matlab*)mfile->data;
 68:   mxArray            *mat;

 71:   if (!ml->rank) {
 72:     PetscInfo1(mfile,"Putting MATLAB array %s\n",name);
 73: #if !defined(PETSC_USE_COMPLEX)
 74:     mat  = mxCreateDoubleMatrix(m,n,mxREAL);
 75: #else
 76:     mat  = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
 77: #endif
 78:     PetscMemcpy(mxGetPr(mat),array,m*n*sizeof(PetscScalar));
 79:     matPutVariable(ml->ep,name,mat);

 81:     PetscInfo1(mfile,"Put MATLAB array %s\n",name);
 82:   }
 83:   return(0);
 84: }

 88: PetscErrorCode  PetscViewerMatlabPutVariable(PetscViewer viewer,const char *name,void *mat)
 89: {
 90:   PetscViewer_Matlab *ml = (PetscViewer_Matlab*)viewer->data;

 93:   matPutVariable(ml->ep,name,(mxArray*)mat);
 94:   return(0);
 95: }

 99: /*@C
100:     PetscViewerMatlabGetArray - Gets a variable from a MATLAB viewer into an array

102:     Not Collective; only processor zero reads in the array

104:     Input Parameters:
105: +    mfile - the MATLAB file viewer
106: .    m,n - the dimensions of the array
107: .    array - the array (represented in one dimension)
108: -    name - the name of the array

110:    Level: advanced

112:      Notes: Only reads in array values on processor 0.

114: @*/
115: PetscErrorCode  PetscViewerMatlabGetArray(PetscViewer mfile,int m,int n,PetscScalar *array,const char *name)
116: {
117:   PetscErrorCode     ierr;
118:   PetscViewer_Matlab *ml = (PetscViewer_Matlab*)mfile->data;
119:   mxArray            *mat;

122:   if (!ml->rank) {
123:     PetscInfo1(mfile,"Getting MATLAB array %s\n",name);
124:     mat  = matGetVariable(ml->ep,name);
125:     if (!mat) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unable to get array %s from matlab",name);
126:     PetscMemcpy(array,mxGetPr(mat),m*n*sizeof(PetscScalar));
127:     PetscInfo1(mfile,"Got MATLAB array %s\n",name);
128:   }
129:   return(0);
130: }

134: PetscErrorCode  PetscViewerFileSetMode_Matlab(PetscViewer viewer,PetscFileMode type)
135: {
136:   PetscViewer_Matlab *vmatlab = (PetscViewer_Matlab*)viewer->data;

139:   vmatlab->btype = type;
140:   return(0);
141: }

143: /*
144:         Actually opens the file
145: */
148: PetscErrorCode  PetscViewerFileSetName_Matlab(PetscViewer viewer,const char name[])
149: {
150:   PetscViewer_Matlab *vmatlab = (PetscViewer_Matlab*)viewer->data;
151:   PetscFileMode      type     = vmatlab->btype;

154:   if (type == (PetscFileMode) -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
155:   if (vmatlab->ep) matClose(vmatlab->ep);

157:   /* only first processor opens file */
158:   if (!vmatlab->rank) {
159:     if (type == FILE_MODE_READ) vmatlab->ep = matOpen(name,"r");
160:     else if (type == FILE_MODE_WRITE || type == FILE_MODE_WRITE) vmatlab->ep = matOpen(name,"w");
161:     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown file type");
162:   }
163:   return(0);
164: }

168: PetscErrorCode PetscViewerDestroy_Matlab(PetscViewer v)
169: {
170:   PetscErrorCode     ierr;
171:   PetscViewer_Matlab *vf = (PetscViewer_Matlab*)v->data;

174:   if (vf->ep) matClose(vf->ep);
175:   PetscFree(vf);
176:   return(0);
177: }

181: PETSC_EXTERN PetscErrorCode PetscViewerCreate_Matlab(PetscViewer viewer)
182: {
183:   PetscErrorCode     ierr;
184:   PetscViewer_Matlab *e;

187:   PetscNewLog(viewer,&e);
188:   MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&e->rank);
189:   e->btype     = (PetscFileMode)-1;
190:   viewer->data = (void*) e;

192:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",PetscViewerFileSetName_Matlab);
193:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_Matlab);

195:   viewer->ops->destroy = PetscViewerDestroy_Matlab;
196:   return(0);
197: }

201: /*@C
202:    PetscViewerMatlabOpen - Opens a Matlab .mat file for output

204:    Collective on MPI_Comm

206:    Input Parameters:
207: +  comm - MPI communicator
208: .  name - name of file
209: -  type - type of file
210: $    FILE_MODE_WRITE - create new file for MATLAB output
211: $    FILE_MODE_READ - open existing file for MATLAB input
212: $    FILE_MODE_WRITE - open existing file for MATLAB output

214:    Output Parameter:
215: .  binv - PetscViewer for MATLAB output to use with the specified file

217:    Level: beginner

219:    Note: This PetscViewer should be destroyed with PetscViewerDestroy().

221:     For writing files it only opens the file on processor 0 in the communicator.

223:      This only saves Vecs it cannot be used to save Mats. We recommend using the PETSCVIEWERBINARY to save objects to be loaded into MATLAB
224:      instead of this routine.

226:    Concepts: MATLAB .mat files
227:    Concepts: PetscViewerMatlab^creating

229: .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(), PETSCVIEWERBINARY, PetscViewerBinaryOpen()
230:           VecView(), MatView(), VecLoad(), MatLoad()
231: @*/
232: PetscErrorCode  PetscViewerMatlabOpen(MPI_Comm comm,const char name[],PetscFileMode type,PetscViewer *binv)
233: {

237:   PetscViewerCreate(comm,binv);
238:   PetscViewerSetType(*binv,PETSCVIEWERMATLAB);
239:   PetscViewerFileSetMode(*binv,type);
240:   PetscViewerFileSetName(*binv,name);
241:   return(0);
242: }

244: static PetscMPIInt Petsc_Viewer_Matlab_keyval = MPI_KEYVAL_INVALID;

248: /*@C
249:      PETSC_VIEWER_MATLAB_ - Creates a Matlab PetscViewer shared by all processors
250:                      in a communicator.

252:      Collective on MPI_Comm

254:      Input Parameter:
255: .    comm - the MPI communicator to share the Matlab PetscViewer

257:      Level: intermediate

259:    Options Database Keys:
260: .    -viewer_matlab_filename <name>

262:    Environmental variables:
263: .   PETSC_VIEWER_MATLAB_FILENAME

265:      Notes:
266:      Unlike almost all other PETSc routines, PETSC_VIEWER_MATLAB_ does not return
267:      an error code.  The matlab PetscViewer is usually used in the form
268: $       XXXView(XXX object,PETSC_VIEWER_MATLAB_(comm));

270:      Use PETSC_VIEWER_SOCKET_() or PetscViewerSocketOpen() to communicator with an interactive MATLAB session.

272: .seealso: PETSC_VIEWER_MATLAB_WORLD, PETSC_VIEWER_MATLAB_SELF, PetscViewerMatlabOpen(), PetscViewerCreate(),
273:           PetscViewerDestroy()
274: @*/
275: PetscViewer  PETSC_VIEWER_MATLAB_(MPI_Comm comm)
276: {
278:   PetscBool      flg;
279:   PetscViewer    viewer;
280:   char           fname[PETSC_MAX_PATH_LEN];
281:   MPI_Comm       ncomm;

284:   PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
285:   if (Petsc_Viewer_Matlab_keyval == MPI_KEYVAL_INVALID) {
286:     MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_Matlab_keyval,0);
287:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
288:   }
289:   MPI_Attr_get(ncomm,Petsc_Viewer_Matlab_keyval,(void**)&viewer,(int*)&flg);
290:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
291:   if (!flg) { /* PetscViewer not yet created */
292:     PetscOptionsGetenv(ncomm,"PETSC_VIEWER_MATLAB_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
293:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
294:     if (!flg) {
295:       PetscStrcpy(fname,"matlaboutput.mat");
296:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
297:     }
298:     PetscViewerMatlabOpen(ncomm,fname,FILE_MODE_WRITE,&viewer);
299:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
300:     PetscObjectRegisterDestroy((PetscObject)viewer);
301:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
302:     MPI_Attr_put(ncomm,Petsc_Viewer_Matlab_keyval,(void*)viewer);
303:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
304:   }
305:   PetscCommDestroy(&ncomm);
306:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
307:   PetscFunctionReturn(viewer);
308: }