Actual source code: vtkv.c

  1: #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>

  3: /*MC
  4:     PetscViewerVTKWriteFunction - functional form used to provide a writer to the `PETSCVIEWERVTK`

  6:      Synopsis:
  7: #include <petscviewer.h>
  8:      PetscViewerVTKWriteFunction(PetscObject object,PetscViewer viewer)

 10:      Input Parameters:
 11: +      object - the PETSc object to be written
 12: -      viewer - viewer it is to be written to

 14:    Level: developer

 16: .seealso: [](sec_viewers), `PETSCVIEWERVTK`, `PetscViewerVTKAddField()`
 17: M*/

 19: /*@C
 20:   PetscViewerVTKAddField - Add a field to the viewer

 22:   Collective

 24:   Input Parameters:
 25: + viewer                      - `PETSCVIEWERVTK`
 26: . dm                          - `DM` on which `Vec` lives
 27: . PetscViewerVTKWriteFunction - function to write this `Vec`
 28: . fieldnum                    - which field of the `DM` to write (`PETSC_DEFAULT` if the whole vector should be written)
 29: . fieldtype                   - Either `PETSC_VTK_POINT_FIELD` or `PETSC_VTK_CELL_FIELD`
 30: . checkdm                     - whether to check for identical dm arguments as fields are added
 31: - vec                         - `Vec` from which to write

 33:   Level: developer

 35:   Note:
 36:   This routine keeps exclusive ownership of the `Vec`. The caller should not use or destroy the `Vec` after calling it.

 38: .seealso: [](sec_viewers), `PETSCVIEWERVTK`, `PetscViewerVTKOpen()`, `DMDAVTKWriteAll()`, `PetscViewerVTKWriteFunction`, `PetscViewerVTKGetDM()`
 39: @*/
 40: PetscErrorCode PetscViewerVTKAddField(PetscViewer viewer, PetscObject dm, PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject, PetscViewer), PetscInt fieldnum, PetscViewerVTKFieldType fieldtype, PetscBool checkdm, PetscObject vec)
 41: {
 42:   PetscFunctionBegin;
 46:   PetscUseMethod(viewer, "PetscViewerVTKAddField_C", (PetscViewer, PetscObject, PetscErrorCode(*)(PetscObject, PetscViewer), PetscInt, PetscViewerVTKFieldType, PetscBool, PetscObject), (viewer, dm, PetscViewerVTKWriteFunction, fieldnum, fieldtype, checkdm, vec));
 47:   PetscFunctionReturn(PETSC_SUCCESS);
 48: }

 50: /*@C
 51:   PetscViewerVTKGetDM - get the `DM` associated with the `PETSCVIEWERVTK` viewer

 53:   Collective

 55:   Input Parameters:
 56: + viewer - `PETSCVIEWERVTK` viewer
 57: - dm     - `DM` associated with the viewer (as a `PetscObject`)

 59:   Level: developer

 61: .seealso: [](sec_viewers), `PETSCVIEWERVTK`, `PetscViewerVTKOpen()`, `DMDAVTKWriteAll()`, `PetscViewerVTKWriteFunction`, `PetscViewerVTKAddField()`
 62: @*/
 63: PetscErrorCode PetscViewerVTKGetDM(PetscViewer viewer, PetscObject *dm)
 64: {
 65:   PetscFunctionBegin;
 67:   PetscUseMethod(viewer, "PetscViewerVTKGetDM_C", (PetscViewer, PetscObject *), (viewer, dm));
 68:   PetscFunctionReturn(PETSC_SUCCESS);
 69: }

 71: static PetscErrorCode PetscViewerDestroy_VTK(PetscViewer viewer)
 72: {
 73:   PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;

 75:   PetscFunctionBegin;
 76:   PetscCall(PetscFree(vtk->filename));
 77:   PetscCall(PetscFree(vtk));
 78:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
 79:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
 80:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
 81:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
 82:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerVTKAddField_C", NULL));
 83:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerVTKGetDM_C", NULL));
 84:   PetscFunctionReturn(PETSC_SUCCESS);
 85: }

 87: static PetscErrorCode PetscViewerFlush_VTK(PetscViewer viewer)
 88: {
 89:   PetscViewer_VTK         *vtk = (PetscViewer_VTK *)viewer->data;
 90:   PetscViewerVTKObjectLink link, next;

 92:   PetscFunctionBegin;
 93:   PetscCheck(!vtk->link || !(!vtk->dm || !vtk->write), PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "No fields or no grid");
 94:   if (vtk->write) PetscCall((*vtk->write)(vtk->dm, viewer));
 95:   for (link = vtk->link; link; link = next) {
 96:     next = link->next;
 97:     PetscCall(PetscObjectDestroy(&link->vec));
 98:     PetscCall(PetscFree(link));
 99:   }
100:   PetscCall(PetscObjectDestroy(&vtk->dm));
101:   vtk->write = NULL;
102:   vtk->link  = NULL;
103:   PetscFunctionReturn(PETSC_SUCCESS);
104: }

106: static PetscErrorCode PetscViewerFileSetName_VTK(PetscViewer viewer, const char name[])
107: {
108:   PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
109:   PetscBool        isvtk, isvts, isvtu, isvtr;
110:   size_t           len;

112:   PetscFunctionBegin;
113:   PetscCall(PetscViewerFlush(viewer));
114:   PetscCall(PetscFree(vtk->filename));
115:   PetscCall(PetscStrlen(name, &len));
116:   if (!len) {
117:     isvtk = PETSC_TRUE;
118:   } else {
119:     PetscCall(PetscStrcasecmp(name + len - 4, ".vtk", &isvtk));
120:     PetscCall(PetscStrcasecmp(name + len - 4, ".vts", &isvts));
121:     PetscCall(PetscStrcasecmp(name + len - 4, ".vtu", &isvtu));
122:     PetscCall(PetscStrcasecmp(name + len - 4, ".vtr", &isvtr));
123:   }
124:   if (isvtk) {
125:     if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_ASCII_VTK_DEPRECATED;
126:     PetscCheck(viewer->format == PETSC_VIEWER_ASCII_VTK_DEPRECATED, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use file '%s' with format %s, should have '.vtk' extension", name, PetscViewerFormats[viewer->format]);
127:   } else if (isvts) {
128:     if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTS;
129:     PetscCheck(viewer->format == PETSC_VIEWER_VTK_VTS, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use file '%s' with format %s, should have '.vts' extension", name, PetscViewerFormats[viewer->format]);
130:   } else if (isvtu) {
131:     if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTU;
132:     PetscCheck(viewer->format == PETSC_VIEWER_VTK_VTU, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use file '%s' with format %s, should have '.vtu' extension", name, PetscViewerFormats[viewer->format]);
133:   } else if (isvtr) {
134:     if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTR;
135:     PetscCheck(viewer->format == PETSC_VIEWER_VTK_VTR, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use file '%s' with format %s, should have '.vtr' extension", name, PetscViewerFormats[viewer->format]);
136:   } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_UNKNOWN_TYPE, "File '%s' has unrecognized extension", name);
137:   PetscCall(PetscStrallocpy(len ? name : "stdout", &vtk->filename));
138:   PetscFunctionReturn(PETSC_SUCCESS);
139: }

141: static PetscErrorCode PetscViewerFileGetName_VTK(PetscViewer viewer, const char **name)
142: {
143:   PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
144:   PetscFunctionBegin;
145:   *name = vtk->filename;
146:   PetscFunctionReturn(PETSC_SUCCESS);
147: }

149: static PetscErrorCode PetscViewerFileSetMode_VTK(PetscViewer viewer, PetscFileMode type)
150: {
151:   PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;

153:   PetscFunctionBegin;
154:   vtk->btype = type;
155:   PetscFunctionReturn(PETSC_SUCCESS);
156: }

158: static PetscErrorCode PetscViewerFileGetMode_VTK(PetscViewer viewer, PetscFileMode *type)
159: {
160:   PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;

162:   PetscFunctionBegin;
163:   *type = vtk->btype;
164:   PetscFunctionReturn(PETSC_SUCCESS);
165: }

167: static PetscErrorCode PetscViewerVTKAddField_VTK(PetscViewer viewer, PetscObject dm, PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject, PetscViewer), PetscInt fieldnum, PetscViewerVTKFieldType fieldtype, PetscBool checkdm, PetscObject vec)
168: {
169:   PetscViewer_VTK         *vtk = (PetscViewer_VTK *)viewer->data;
170:   PetscViewerVTKObjectLink link, tail = vtk->link;

172:   PetscFunctionBegin;
173:   if (vtk->dm) {
174:     PetscCheck(!checkdm || dm == vtk->dm, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Refusing to write a field from more than one grid to the same VTK file. Set checkdm = PETSC_FALSE to skip this check.");
175:   } else {
176:     PetscCall(PetscObjectReference(dm));
177:     vtk->dm = dm;
178:   }
179:   vtk->write = PetscViewerVTKWriteFunction;
180:   PetscCall(PetscNew(&link));
181:   link->ft    = fieldtype;
182:   link->vec   = vec;
183:   link->field = fieldnum;
184:   link->next  = NULL;
185:   /* Append to list */
186:   if (tail) {
187:     while (tail->next) tail = tail->next;
188:     tail->next = link;
189:   } else vtk->link = link;
190:   PetscFunctionReturn(PETSC_SUCCESS);
191: }

193: static PetscErrorCode PetscViewerVTKGetDM_VTK(PetscViewer viewer, PetscObject *dm)
194: {
195:   PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;

197:   PetscFunctionBegin;
198:   *dm = vtk->dm;
199:   PetscFunctionReturn(PETSC_SUCCESS);
200: }

202: /*MC
203:    PETSCVIEWERVTK - A viewer that writes to a VTK file

205:   Level: beginner

207: .seealso: [](sec_viewers), `PetscViewerVTKOpen()`, `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
208:           `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
209:           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
210:           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
211: M*/

213: PETSC_EXTERN PetscErrorCode PetscViewerCreate_VTK(PetscViewer v)
214: {
215:   PetscViewer_VTK *vtk;

217:   PetscFunctionBegin;
218:   PetscCall(PetscNew(&vtk));

220:   v->data         = (void *)vtk;
221:   v->ops->destroy = PetscViewerDestroy_VTK;
222:   v->ops->flush   = PetscViewerFlush_VTK;
223:   vtk->btype      = FILE_MODE_UNDEFINED;
224:   vtk->filename   = NULL;

226:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_VTK));
227:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_VTK));
228:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_VTK));
229:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_VTK));
230:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerVTKAddField_C", PetscViewerVTKAddField_VTK));
231:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerVTKGetDM_C", PetscViewerVTKGetDM_VTK));
232:   PetscFunctionReturn(PETSC_SUCCESS);
233: }

235: /*@C
236:   PetscViewerVTKOpen - Opens a `PETSCVIEWERVTK` viewer file.

238:   Collective

240:   Input Parameters:
241: + comm - MPI communicator
242: . name - name of file
243: - type - type of file
244: .vb
245:     FILE_MODE_WRITE - create new file for binary output
246:     FILE_MODE_READ - open existing file for binary input (not currently supported)
247:     FILE_MODE_APPEND - open existing file for binary output (not currently supported)
248: .ve

250:   Output Parameter:
251: . vtk - `PetscViewer` for VTK input/output to use with the specified file

253:   Level: beginner

255: .seealso: [](sec_viewers), `PETSCVIEWERVTK`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
256:           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`,
257:           `PetscFileMode`, `PetscViewer`
258: @*/
259: PetscErrorCode PetscViewerVTKOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *vtk)
260: {
261:   PetscFunctionBegin;
262:   PetscCall(PetscViewerCreate(comm, vtk));
263:   PetscCall(PetscViewerSetType(*vtk, PETSCVIEWERVTK));
264:   PetscCall(PetscViewerFileSetMode(*vtk, type));
265:   PetscCall(PetscViewerFileSetName(*vtk, name));
266:   PetscFunctionReturn(PETSC_SUCCESS);
267: }

269: /*@C
270:   PetscViewerVTKFWrite - write binary data preceded by 32-bit int length (in bytes), does not do byte swapping.

272:   Logically Collective

274:   Input Parameters:
275: + viewer - logically collective viewer, data written from rank 0
276: . fp     - file pointer valid on rank 0
277: . data   - data pointer valid on rank 0
278: . n      - number of data items
279: - dtype  - data type

281:   Level: developer

283:   Note:
284:   If `PetscScalar` is `__float128` then the binary files are written in double precision

286: .seealso: [](sec_viewers), `PETSCVIEWERVTK`, `DMDAVTKWriteAll()`, `DMPlexVTKWriteAll()`, `PetscViewerPushFormat()`, `PetscViewerVTKOpen()`, `PetscBinaryWrite()`
287: @*/
288: PetscErrorCode PetscViewerVTKFWrite(PetscViewer viewer, FILE *fp, const void *data, PetscInt n, MPI_Datatype dtype)
289: {
290:   PetscMPIInt  rank;
291:   MPI_Datatype vdtype = dtype;
292: #if defined(PETSC_USE_REAL___FLOAT128)
293:   double    *tmp;
294:   PetscInt   i;
295:   PetscReal *ttmp = (PetscReal *)data;
296: #endif

298:   PetscFunctionBegin;
299:   PetscCheck(n >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Trying to write a negative amount of data %" PetscInt_FMT, n);
300:   if (!n) PetscFunctionReturn(PETSC_SUCCESS);
301:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
302:   if (rank == 0) {
303:     size_t      count;
304:     PetscMPIInt dsize;
305:     PetscInt64  bytes;

307: #if defined(PETSC_USE_REAL___FLOAT128)
308:     if (dtype == MPIU___FLOAT128) {
309:       PetscCall(PetscMalloc1(n, &tmp));
310:       for (i = 0; i < n; i++) tmp[i] = ttmp[i];
311:       data   = (void *)tmp;
312:       vdtype = MPI_DOUBLE;
313:     }
314: #endif
315:     PetscCallMPI(MPI_Type_size(vdtype, &dsize));
316:     bytes = (PetscInt64)dsize * n;

318:     count = fwrite(&bytes, sizeof(bytes), 1, fp);
319:     PetscCheck(count == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_WRITE, "Error writing byte count");
320:     count = fwrite(data, dsize, (size_t)n, fp);
321:     PetscCheck((PetscInt)count == n, PETSC_COMM_SELF, PETSC_ERR_FILE_WRITE, "Wrote %" PetscInt_FMT "/%" PetscInt_FMT " array members of size %d", (PetscInt)count, n, dsize);
322: #if defined(PETSC_USE_REAL___FLOAT128)
323:     if (dtype == MPIU___FLOAT128) PetscCall(PetscFree(tmp));
324: #endif
325:   }
326:   PetscFunctionReturn(PETSC_SUCCESS);
327: }