Actual source code: vmatlab.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2:  #include <petsc/private/viewerimpl.h>
  3: #include <mat.h>


  6: typedef struct {
  7:   MATFile       *ep;
  8:   PetscMPIInt   rank;
  9:   PetscFileMode btype;
 10: } PetscViewer_Matlab;

 12: /*@C
 13:     PetscViewerMatlabPutArray - Puts an array into the MATLAB viewer.

 15:       Not collective: only processor zero saves the array

 17:     Input Parameters:
 18: +    mfile - the viewer
 19: .    m,n - the dimensions of the array
 20: .    array - the array (represented in one dimension)
 21: -    name - the name of the array

 23:    Level: advanced

 25:      Notes:
 26:     Only writes array values on processor 0.

 28: @*/
 29: PetscErrorCode  PetscViewerMatlabPutArray(PetscViewer mfile,int m,int n,const PetscScalar *array,const char *name)
 30: {
 31:   PetscErrorCode     ierr;
 32:   PetscViewer_Matlab *ml = (PetscViewer_Matlab*)mfile->data;
 33:   mxArray            *mat;

 36:   if (!ml->rank) {
 37:     PetscInfo1(mfile,"Putting MATLAB array %s\n",name);
 38: #if !defined(PETSC_USE_COMPLEX)
 39:     mat  = mxCreateDoubleMatrix(m,n,mxREAL);
 40: #else
 41:     mat  = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
 42: #endif
 43:     PetscMemcpy(mxGetPr(mat),array,m*n*sizeof(PetscScalar));
 44:     matPutVariable(ml->ep,name,mat);

 46:     PetscInfo1(mfile,"Put MATLAB array %s\n",name);
 47:   }
 48:   return(0);
 49: }

 51: PetscErrorCode  PetscViewerMatlabPutVariable(PetscViewer viewer,const char *name,void *mat)
 52: {
 53:   PetscViewer_Matlab *ml = (PetscViewer_Matlab*)viewer->data;

 56:   matPutVariable(ml->ep,name,(mxArray*)mat);
 57:   return(0);
 58: }

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

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

 65:     Input Parameters:
 66: +    mfile - the MATLAB file viewer
 67: .    m,n - the dimensions of the array
 68: .    array - the array (represented in one dimension)
 69: -    name - the name of the array

 71:    Level: advanced

 73:      Notes:
 74:     Only reads in array values on processor 0.

 76: @*/
 77: PetscErrorCode  PetscViewerMatlabGetArray(PetscViewer mfile,int m,int n,PetscScalar *array,const char *name)
 78: {
 79:   PetscErrorCode     ierr;
 80:   PetscViewer_Matlab *ml = (PetscViewer_Matlab*)mfile->data;
 81:   mxArray            *mat;

 84:   if (!ml->rank) {
 85:     PetscInfo1(mfile,"Getting MATLAB array %s\n",name);
 86:     mat  = matGetVariable(ml->ep,name);
 87:     if (!mat) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unable to get array %s from matlab",name);
 88:     PetscMemcpy(array,mxGetPr(mat),m*n*sizeof(PetscScalar));
 89:     PetscInfo1(mfile,"Got MATLAB array %s\n",name);
 90:   }
 91:   return(0);
 92: }

 94: PetscErrorCode  PetscViewerFileSetMode_Matlab(PetscViewer viewer,PetscFileMode type)
 95: {
 96:   PetscViewer_Matlab *vmatlab = (PetscViewer_Matlab*)viewer->data;

 99:   vmatlab->btype = type;
100:   return(0);
101: }

103: /*
104:         Actually opens the file
105: */
106: PetscErrorCode  PetscViewerFileSetName_Matlab(PetscViewer viewer,const char name[])
107: {
108:   PetscViewer_Matlab *vmatlab = (PetscViewer_Matlab*)viewer->data;
109:   PetscFileMode      type     = vmatlab->btype;

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

115:   /* only first processor opens file */
116:   if (!vmatlab->rank) {
117:     if (type == FILE_MODE_READ) vmatlab->ep = matOpen(name,"r");
118:     else if (type == FILE_MODE_WRITE || type == FILE_MODE_WRITE) vmatlab->ep = matOpen(name,"w");
119:     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown file type");
120:   }
121:   return(0);
122: }

124: PetscErrorCode PetscViewerDestroy_Matlab(PetscViewer v)
125: {
126:   PetscErrorCode     ierr;
127:   PetscViewer_Matlab *vf = (PetscViewer_Matlab*)v->data;

130:   if (vf->ep) matClose(vf->ep);
131:   PetscFree(vf);
132:   return(0);
133: }

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

139:    Level: intermediate

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

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

149: $             In your PETSc C/C++ code (assuming a two dimensional DMDA with one degree of freedom per node)
150: $                PetscObjectSetName((PetscObject)x,"x");
151: $                VecView(x,PETSC_VIEWER_MATLAB_WORLD);
152: $                PetscObjectSetName((PetscObject)da,"da");
153: $                DMView(x,PETSC_VIEWER_MATLAB_WORLD);
154: $             Then from MATLAB
155: $                load('matlaboutput.mat')   % matlaboutput.mat is the default filename
156: $                xnew = zeros(da.n,da.m);
157: $                xnew(:) = x;    % reshape one dimensional vector back to two dimensions

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

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

164: .seealso:  PETSC_VIEWER_MATLAB_(),PETSC_VIEWER_MATLAB_SELF, PETSC_VIEWER_MATLAB_WORLD,PetscViewerCreate(),
165:            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERBINARY, PETSCVIEWERASCII, PETSCVIEWERDRAW,
166:            PETSC_VIEWER_STDOUT_(), PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat

168: M*/
169: PETSC_EXTERN PetscErrorCode PetscViewerCreate_Matlab(PetscViewer viewer)
170: {
171:   PetscErrorCode     ierr;
172:   PetscViewer_Matlab *e;

175:   PetscNewLog(viewer,&e);
176:   MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&e->rank);
177:   e->btype     = (PetscFileMode)-1;
178:   viewer->data = (void*) e;

180:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",PetscViewerFileSetName_Matlab);
181:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_Matlab);

183:   viewer->ops->destroy = PetscViewerDestroy_Matlab;
184:   return(0);
185: }

187: /*@C
188:    PetscViewerMatlabOpen - Opens a Matlab .mat file for output

190:    Collective on MPI_Comm

192:    Input Parameters:
193: +  comm - MPI communicator
194: .  name - name of file
195: -  type - type of file
196: $    FILE_MODE_WRITE - create new file for MATLAB output
197: $    FILE_MODE_READ - open existing file for MATLAB input
198: $    FILE_MODE_WRITE - open existing file for MATLAB output

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

203:    Level: beginner

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

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

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

212:    Concepts: MATLAB .mat files
213:    Concepts: PetscViewerMatlab^creating

215: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PETSCVIEWERBINARY, PetscViewerBinaryOpen()
216:           VecView(), MatView(), VecLoad(), MatLoad()
217: @*/
218: PetscErrorCode  PetscViewerMatlabOpen(MPI_Comm comm,const char name[],PetscFileMode type,PetscViewer *binv)
219: {

223:   PetscViewerCreate(comm,binv);
224:   PetscViewerSetType(*binv,PETSCVIEWERMATLAB);
225:   PetscViewerFileSetMode(*binv,type);
226:   PetscViewerFileSetName(*binv,name);
227:   return(0);
228: }

230: static PetscMPIInt Petsc_Viewer_Matlab_keyval = MPI_KEYVAL_INVALID;

232: /*@C
233:      PETSC_VIEWER_MATLAB_ - Creates a Matlab PetscViewer shared by all processors
234:                      in a communicator.

236:      Collective on MPI_Comm

238:      Input Parameter:
239: .    comm - the MPI communicator to share the Matlab PetscViewer

241:      Level: intermediate

243:    Options Database Keys:
244: .    -viewer_matlab_filename <name>

246:    Environmental variables:
247: .   PETSC_VIEWER_MATLAB_FILENAME

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

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

256: .seealso: PETSC_VIEWER_MATLAB_WORLD, PETSC_VIEWER_MATLAB_SELF, PetscViewerMatlabOpen(), PetscViewerCreate(),
257:           PetscViewerDestroy()
258: @*/
259: PetscViewer  PETSC_VIEWER_MATLAB_(MPI_Comm comm)
260: {
262:   PetscBool      flg;
263:   PetscViewer    viewer;
264:   char           fname[PETSC_MAX_PATH_LEN];
265:   MPI_Comm       ncomm;

268:   PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
269:   if (Petsc_Viewer_Matlab_keyval == MPI_KEYVAL_INVALID) {
270:     MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_Matlab_keyval,0);
271:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
272:   }
273:   MPI_Comm_get_attr(ncomm,Petsc_Viewer_Matlab_keyval,(void**)&viewer,(int*)&flg);
274:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
275:   if (!flg) { /* PetscViewer not yet created */
276:     PetscOptionsGetenv(ncomm,"PETSC_VIEWER_MATLAB_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
277:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
278:     if (!flg) {
279:       PetscStrcpy(fname,"matlaboutput.mat");
280:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
281:     }
282:     PetscViewerMatlabOpen(ncomm,fname,FILE_MODE_WRITE,&viewer);
283:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
284:     PetscObjectRegisterDestroy((PetscObject)viewer);
285:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
286:     MPI_Comm_set_attr(ncomm,Petsc_Viewer_Matlab_keyval,(void*)viewer);
287:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
288:   }
289:   PetscCommDestroy(&ncomm);
290:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_MATLAB_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
291:   PetscFunctionReturn(viewer);
292: }