Actual source code: plexvtk.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: #define PETSCDM_DLL
  2:  #include <petsc/private/dmpleximpl.h>
  3:  #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>

  5: PetscErrorCode DMPlexVTKGetCellType_Internal(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType)
  6: {
  8:   *cellType = -1;
  9:   switch (dim) {
 10:   case 0:
 11:     switch (corners) {
 12:     case 1:
 13:       *cellType = 1; /* VTK_VERTEX */
 14:       break;
 15:     default:
 16:       break;
 17:     }
 18:     break;
 19:   case 1:
 20:     switch (corners) {
 21:     case 2:
 22:       *cellType = 3; /* VTK_LINE */
 23:       break;
 24:     case 3:
 25:       *cellType = 21; /* VTK_QUADRATIC_EDGE */
 26:       break;
 27:     default:
 28:       break;
 29:     }
 30:     break;
 31:   case 2:
 32:     switch (corners) {
 33:     case 3:
 34:       *cellType = 5; /* VTK_TRIANGLE */
 35:       break;
 36:     case 4:
 37:       *cellType = 9; /* VTK_QUAD */
 38:       break;
 39:     case 6:
 40:       *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */
 41:       break;
 42:     case 9:
 43:       *cellType = 23; /* VTK_QUADRATIC_QUAD */
 44:       break;
 45:     default:
 46:       break;
 47:     }
 48:     break;
 49:   case 3:
 50:     switch (corners) {
 51:     case 4:
 52:       *cellType = 10; /* VTK_TETRA */
 53:       break;
 54:     case 6:
 55:       *cellType = 13; /* VTK_WEDGE */
 56:       break;
 57:     case 8:
 58:       *cellType = 12; /* VTK_HEXAHEDRON */
 59:       break;
 60:     case 10:
 61:       *cellType = 24; /* VTK_QUADRATIC_TETRA */
 62:       break;
 63:     case 27:
 64:       *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */
 65:       break;
 66:     default:
 67:       break;
 68:     }
 69:   }
 70:   return(0);
 71: }

 73: static PetscErrorCode DMPlexVTKWriteCells_ASCII(DM dm, FILE *fp, PetscInt *totalCells)
 74: {
 75:   MPI_Comm       comm;
 76:   DMLabel        label;
 77:   IS             globalVertexNumbers = NULL;
 78:   const PetscInt *gvertex;
 79:   PetscInt       dim;
 80:   PetscInt       numCorners = 0, totCorners = 0, maxCorners, *corners;
 81:   PetscInt       numCells   = 0, totCells   = 0, maxCells, cellHeight;
 82:   PetscInt       numLabelCells, maxLabelCells, cStart, cEnd, c, vStart, vEnd, v;
 83:   PetscMPIInt    size, rank, proc, tag;

 87:   PetscObjectGetComm((PetscObject)dm,&comm);
 88:   PetscCommGetNewTag(comm, &tag);
 89:   MPI_Comm_size(comm, &size);
 90:   MPI_Comm_rank(comm, &rank);
 91:   DMGetDimension(dm, &dim);
 92:   DMPlexGetVTKCellHeight(dm, &cellHeight);
 93:   DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
 94:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
 95:   DMGetLabel(dm, "vtk", &label);
 96:   DMGetStratumSize(dm, "vtk", 1, &numLabelCells);
 97:   MPIU_Allreduce(&numLabelCells, &maxLabelCells, 1, MPIU_INT, MPI_MAX, comm);
 98:   if (!maxLabelCells) label = NULL;
 99:   for (c = cStart; c < cEnd; ++c) {
100:     PetscInt *closure = NULL;
101:     PetscInt closureSize, value;

103:     if (label) {
104:       DMLabelGetValue(label, c, &value);
105:       if (value != 1) continue;
106:     }
107:     DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
108:     for (v = 0; v < closureSize*2; v += 2) {
109:       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners;
110:     }
111:     DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
112:     ++numCells;
113:   }
114:   maxCells = numCells;
115:   MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);
116:   MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);
117:   MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);
118:   MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);
119:   DMPlexGetVertexNumbering(dm, &globalVertexNumbers);
120:   ISGetIndices(globalVertexNumbers, &gvertex);
121:   PetscMalloc1(maxCells, &corners);
122:   PetscFPrintf(comm, fp, "CELLS %D %D\n", totCells, totCorners+totCells);
123:   if (!rank) {
124:     PetscInt *remoteVertices;
125:     PetscInt *vertices;

127:     PetscMalloc1(maxCorners, &vertices);
128:     for (c = cStart, numCells = 0; c < cEnd; ++c) {
129:       PetscInt *closure = NULL;
130:       PetscInt closureSize, value, nC = 0;

132:       if (label) {
133:         DMLabelGetValue(label, c, &value);
134:         if (value != 1) continue;
135:       }
136:       DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
137:       for (v = 0; v < closureSize*2; v += 2) {
138:         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
139:           const PetscInt gv = gvertex[closure[v] - vStart];
140:           vertices[nC++] = gv < 0 ? -(gv+1) : gv;
141:         }
142:       }
143:       DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
144:       DMPlexReorderCell(dm, c, vertices);
145:       corners[numCells++] = nC;
146:       PetscFPrintf(comm, fp, "%D ", nC);
147:       for (v = 0; v < nC; ++v) {
148:         PetscFPrintf(comm, fp, " %D", vertices[v]);
149:       }
150:       PetscFPrintf(comm, fp, "\n");
151:     }
152:     if (size > 1) {PetscMalloc1(maxCorners+maxCells, &remoteVertices);}
153:     for (proc = 1; proc < size; ++proc) {
154:       MPI_Status status;

156:       MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);
157:       MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);
158:       for (c = 0; c < numCorners;) {
159:         PetscInt nC = remoteVertices[c++];

161:         for (v = 0; v < nC; ++v, ++c) {
162:           vertices[v] = remoteVertices[c];
163:         }
164:         PetscFPrintf(comm, fp, "%D ", nC);
165:         for (v = 0; v < nC; ++v) {
166:           PetscFPrintf(comm, fp, " %D", vertices[v]);
167:         }
168:         PetscFPrintf(comm, fp, "\n");
169:       }
170:     }
171:     if (size > 1) {PetscFree(remoteVertices);}
172:     PetscFree(vertices);
173:   } else {
174:     PetscInt *localVertices, numSend = numCells+numCorners, k = 0;

176:     PetscMalloc1(numSend, &localVertices);
177:     for (c = cStart, numCells = 0; c < cEnd; ++c) {
178:       PetscInt *closure = NULL;
179:       PetscInt closureSize, value, nC = 0;

181:       if (label) {
182:         DMLabelGetValue(label, c, &value);
183:         if (value != 1) continue;
184:       }
185:       DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
186:       for (v = 0; v < closureSize*2; v += 2) {
187:         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
188:           const PetscInt gv = gvertex[closure[v] - vStart];
189:           closure[nC++] = gv < 0 ? -(gv+1) : gv;
190:         }
191:       }
192:       corners[numCells++] = nC;
193:       localVertices[k++]  = nC;
194:       for (v = 0; v < nC; ++v, ++k) {
195:         localVertices[k] = closure[v];
196:       }
197:       DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
198:       DMPlexReorderCell(dm, c, localVertices+k-nC);
199:     }
200:     if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %D should be %D", k, numSend);
201:     MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);
202:     MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);
203:     PetscFree(localVertices);
204:   }
205:   ISRestoreIndices(globalVertexNumbers, &gvertex);
206:   PetscFPrintf(comm, fp, "CELL_TYPES %D\n", totCells);
207:   if (!rank) {
208:     PetscInt cellType;

210:     for (c = 0; c < numCells; ++c) {
211:       DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType);
212:       PetscFPrintf(comm, fp, "%D\n", cellType);
213:     }
214:     for (proc = 1; proc < size; ++proc) {
215:       MPI_Status status;

217:       MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);
218:       MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);
219:       for (c = 0; c < numCells; ++c) {
220:         DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType);
221:         PetscFPrintf(comm, fp, "%D\n", cellType);
222:       }
223:     }
224:   } else {
225:     MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);
226:     MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);
227:   }
228:   PetscFree(corners);
229:   *totalCells = totCells;
230:   return(0);
231: }

233: static PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp)
234: {
235:   MPI_Comm       comm;
236:   PetscInt       numCells = 0, cellHeight;
237:   PetscInt       numLabelCells, cStart, cEnd, c;
238:   PetscMPIInt    size, rank, proc, tag;
239:   PetscBool      hasLabel;

243:   PetscObjectGetComm((PetscObject)dm,&comm);
244:   PetscCommGetNewTag(comm, &tag);
245:   MPI_Comm_size(comm, &size);
246:   MPI_Comm_rank(comm, &rank);
247:   DMPlexGetVTKCellHeight(dm, &cellHeight);
248:   DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
249:   DMGetStratumSize(dm, "vtk", 1, &numLabelCells);
250:   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
251:   for (c = cStart; c < cEnd; ++c) {
252:     if (hasLabel) {
253:       PetscInt value;

255:       DMGetLabelValue(dm, "vtk", c, &value);
256:       if (value != 1) continue;
257:     }
258:     ++numCells;
259:   }
260:   if (!rank) {
261:     for (c = 0; c < numCells; ++c) {PetscFPrintf(comm, fp, "%d\n", rank);}
262:     for (proc = 1; proc < size; ++proc) {
263:       MPI_Status status;

265:       MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);
266:       for (c = 0; c < numCells; ++c) {PetscFPrintf(comm, fp, "%d\n", proc);}
267:     }
268:   } else {
269:     MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);
270:   }
271:   return(0);
272: }

274: #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
275: typedef double PetscVTKReal;
276: #elif defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
277: typedef float PetscVTKReal;
278: #else
279: typedef PetscReal PetscVTKReal;
280: #endif

282: static PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale, PetscInt imag)
283: {
284:   MPI_Comm           comm;
285:   const MPI_Datatype mpiType = MPIU_SCALAR;
286:   PetscScalar        *array;
287:   PetscInt           numDof = 0, maxDof;
288:   PetscInt           numLabelCells, cellHeight, cStart, cEnd, numLabelVertices, vStart, vEnd, pStart, pEnd, p;
289:   PetscMPIInt        size, rank, proc, tag;
290:   PetscBool          hasLabel;
291:   PetscErrorCode     ierr;

294:   PetscObjectGetComm((PetscObject)dm,&comm);
297:   if (precision < 0) precision = 6;
298:   PetscCommGetNewTag(comm, &tag);
299:   MPI_Comm_size(comm, &size);
300:   MPI_Comm_rank(comm, &rank);
301:   PetscSectionGetChart(section, &pStart, &pEnd);
302:   /* VTK only wants the values at cells or vertices */
303:   DMPlexGetVTKCellHeight(dm, &cellHeight);
304:   DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
305:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
306:   pStart   = PetscMax(PetscMin(cStart, vStart), pStart);
307:   pEnd     = PetscMin(PetscMax(cEnd,   vEnd),   pEnd);
308:   DMGetStratumSize(dm, "vtk", 1, &numLabelCells);
309:   DMGetStratumSize(dm, "vtk", 2, &numLabelVertices);
310:   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
311:   for (p = pStart; p < pEnd; ++p) {
312:     /* Reject points not either cells or vertices */
313:     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
314:     if (hasLabel) {
315:       PetscInt value;

317:       if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
318:           ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
319:         DMGetLabelValue(dm, "vtk", p, &value);
320:         if (value != 1) continue;
321:       }
322:     }
323:     PetscSectionGetDof(section, p, &numDof);
324:     if (numDof) break;
325:   }
326:   MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);
327:   enforceDof = PetscMax(enforceDof, maxDof);
328:   VecGetArray(v, &array);
329:   if (!rank) {
330:     PetscVTKReal dval;
331:     PetscScalar  val;
332:     char formatString[8];

334:     PetscSNPrintf(formatString, 8, "%%.%de", precision);
335:     for (p = pStart; p < pEnd; ++p) {
336:       /* Here we lose a way to filter points by keeping them out of the Numbering */
337:       PetscInt dof, off, goff, d;

339:       /* Reject points not either cells or vertices */
340:       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
341:       if (hasLabel) {
342:         PetscInt value;

344:         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
345:             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
346:           DMGetLabelValue(dm, "vtk", p, &value);
347:           if (value != 1) continue;
348:         }
349:       }
350:       PetscSectionGetDof(section, p, &dof);
351:       PetscSectionGetOffset(section, p, &off);
352:       PetscSectionGetOffset(globalSection, p, &goff);
353:       if (dof && goff >= 0) {
354:         for (d = 0; d < dof; d++) {
355:           if (d > 0) {
356:             PetscFPrintf(comm, fp, " ");
357:           }
358:           val = array[off+d];
359:           dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
360:           PetscFPrintf(comm, fp, formatString, dval);
361:         }
362:         for (d = dof; d < enforceDof; d++) {
363:           PetscFPrintf(comm, fp, " 0.0");
364:         }
365:         PetscFPrintf(comm, fp, "\n");
366:       }
367:     }
368:     for (proc = 1; proc < size; ++proc) {
369:       PetscScalar *remoteValues;
370:       PetscInt    size = 0, d;
371:       MPI_Status  status;

373:       MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);
374:       PetscMalloc1(size, &remoteValues);
375:       MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);
376:       for (p = 0; p < size/maxDof; ++p) {
377:         for (d = 0; d < maxDof; ++d) {
378:           if (d > 0) {
379:             PetscFPrintf(comm, fp, " ");
380:           }
381:           val = remoteValues[p*maxDof+d];
382:           dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
383:           PetscFPrintf(comm, fp, formatString, dval);
384:         }
385:         for (d = maxDof; d < enforceDof; ++d) {
386:           PetscFPrintf(comm, fp, " 0.0");
387:         }
388:         PetscFPrintf(comm, fp, "\n");
389:       }
390:       PetscFree(remoteValues);
391:     }
392:   } else {
393:     PetscScalar *localValues;
394:     PetscInt    size, k = 0;

396:     PetscSectionGetStorageSize(section, &size);
397:     PetscMalloc1(size, &localValues);
398:     for (p = pStart; p < pEnd; ++p) {
399:       PetscInt dof, off, goff, d;

401:       /* Reject points not either cells or vertices */
402:       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
403:       if (hasLabel) {
404:         PetscInt value;

406:         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
407:             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
408:           DMGetLabelValue(dm, "vtk", p, &value);
409:           if (value != 1) continue;
410:         }
411:       }
412:       PetscSectionGetDof(section, p, &dof);
413:       PetscSectionGetOffset(section, p, &off);
414:       PetscSectionGetOffset(globalSection, p, &goff);
415:       if (goff >= 0) {
416:         for (d = 0; d < dof; ++d) {
417:           localValues[k++] = array[off+d];
418:         }
419:       }
420:     }
421:     MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);
422:     MPI_Send(localValues, k, mpiType, 0, tag, comm);
423:     PetscFree(localValues);
424:   }
425:   VecRestoreArray(v, &array);
426:   return(0);
427: }

429: static PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale, PetscBool nameComplex, PetscInt imag)
430: {
431:   MPI_Comm       comm;
432:   PetscInt       numDof = 0, maxDof;
433:   PetscInt       pStart, pEnd, p;

437:   PetscObjectGetComm((PetscObject)dm,&comm);
438:   PetscSectionGetChart(section, &pStart, &pEnd);
439:   for (p = pStart; p < pEnd; ++p) {
440:     PetscSectionGetDof(section, p, &numDof);
441:     if (numDof) break;
442:   }
443:   numDof = PetscMax(numDof, enforceDof);
444:   MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));
445:   if (!name) name = "Unknown";
446:   if (maxDof == 3) {
447:     if (nameComplex) {
448:       PetscFPrintf(comm, fp, "VECTORS %s.%s double\n", name, imag ? "Im" : "Re");
449:     } else {
450:       PetscFPrintf(comm, fp, "VECTORS %s double\n", name);
451:     }
452:   } else {
453:     if (nameComplex) {
454:       PetscFPrintf(comm, fp, "SCALARS %s.%s double %D\n", name, imag ? "Im" : "Re", maxDof);
455:     } else {
456:       PetscFPrintf(comm, fp, "SCALARS %s double %D\n", name, maxDof);
457:     }
458:     PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");
459:   }
460:   DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale, imag);
461:   return(0);
462: }

464: static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
465: {
466:   MPI_Comm                 comm;
467:   PetscViewer_VTK          *vtk = (PetscViewer_VTK*) viewer->data;
468:   FILE                     *fp;
469:   PetscViewerVTKObjectLink link;
470:   PetscSection             coordSection, globalCoordSection;
471:   PetscLayout              vLayout;
472:   Vec                      coordinates;
473:   PetscReal                lengthScale;
474:   PetscInt                 totVertices, totCells = 0, loops_per_scalar, l;
475:   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE, localized, writeComplex;
476:   PetscErrorCode           ierr;

479: #if defined(PETSC_USE_COMPLEX)
480:   loops_per_scalar = 2;
481:   writeComplex = PETSC_TRUE;
482: #else
483:   loops_per_scalar = 1;
484:   writeComplex = PETSC_FALSE;
485: #endif
486:   DMGetCoordinatesLocalized(dm,&localized);
487:   if (localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"VTK output with localized coordinates not yet supported");
488:   PetscObjectGetComm((PetscObject)dm,&comm);
489:   PetscFOpen(comm, vtk->filename, "wb", &fp);
490:   PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");
491:   PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");
492:   PetscFPrintf(comm, fp, "ASCII\n");
493:   PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");
494:   /* Vertices */
495:   DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);
496:   DMGetCoordinateSection(dm, &coordSection);
497:   PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalCoordSection);
498:   DMGetCoordinatesLocal(dm, &coordinates);
499:   PetscSectionGetPointLayout(PetscObjectComm((PetscObject)dm), globalCoordSection, &vLayout);
500:   PetscLayoutGetSize(vLayout, &totVertices);
501:   PetscFPrintf(comm, fp, "POINTS %D double\n", totVertices);
502:   DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale, 0);
503:   /* Cells */
504:   DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);
505:   /* Vertex fields */
506:   for (link = vtk->link; link; link = link->next) {
507:     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
508:     if ((link->ft == PETSC_VTK_CELL_FIELD)  || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD))  hasCell  = PETSC_TRUE;
509:   }
510:   if (hasPoint) {
511:     PetscFPrintf(comm, fp, "POINT_DATA %D\n", totVertices);
512:     for (link = vtk->link; link; link = link->next) {
513:       Vec          X = (Vec) link->vec;
514:       PetscSection section = NULL, globalSection, newSection = NULL;
515:       char         namebuf[256];
516:       const char   *name;
517:       PetscInt     enforceDof = PETSC_DETERMINE;

519:       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
520:       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
521:       PetscObjectGetName(link->vec, &name);
522:       PetscObjectQuery(link->vec, "section", (PetscObject*) &section);
523:       if (!section) {
524:         DM           dmX;

526:         VecGetDM(X, &dmX);
527:         if (dmX) {
528:           DMLabel  subpointMap, subpointMapX;
529:           PetscInt dim, dimX, pStart, pEnd, qStart, qEnd;

531:           DMGetLocalSection(dmX, &section);
532:           /* Here is where we check whether dmX is a submesh of dm */
533:           DMGetDimension(dm,  &dim);
534:           DMGetDimension(dmX, &dimX);
535:           DMPlexGetChart(dm,  &pStart, &pEnd);
536:           DMPlexGetChart(dmX, &qStart, &qEnd);
537:           DMPlexGetSubpointMap(dm,  &subpointMap);
538:           DMPlexGetSubpointMap(dmX, &subpointMapX);
539:           if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) {
540:             const PetscInt *ind = NULL;
541:             IS              subpointIS;
542:             PetscInt        n = 0, q;

544:             PetscSectionGetChart(section, &qStart, &qEnd);
545:             DMPlexCreateSubpointIS(dm, &subpointIS);
546:             if (subpointIS) {
547:               ISGetLocalSize(subpointIS, &n);
548:               ISGetIndices(subpointIS, &ind);
549:             }
550:             PetscSectionCreate(comm, &newSection);
551:             PetscSectionSetChart(newSection, pStart, pEnd);
552:             for (q = qStart; q < qEnd; ++q) {
553:               PetscInt dof, off, p;

555:               PetscSectionGetDof(section, q, &dof);
556:               if (dof) {
557:                 PetscFindInt(q, n, ind, &p);
558:                 if (p >= pStart) {
559:                   PetscSectionSetDof(newSection, p, dof);
560:                   PetscSectionGetOffset(section, q, &off);
561:                   PetscSectionSetOffset(newSection, p, off);
562:                 }
563:               }
564:             }
565:             if (subpointIS) {
566:               ISRestoreIndices(subpointIS, &ind);
567:               ISDestroy(&subpointIS);
568:             }
569:             /* No need to setup section */
570:             section = newSection;
571:           }
572:         }
573:       }
574:       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name);
575:       if (link->field >= 0) {
576:         const char *fieldname;

578:         PetscSectionGetFieldName(section, link->field, &fieldname);
579:         PetscSectionGetField(section, link->field, &section);
580:         if (fieldname) {
581:           PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname);
582:         } else {
583:           PetscSNPrintf(namebuf, sizeof(namebuf), "%s%D", name, link->field);
584:         }
585:       } else {
586:         PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name);
587:       }
588:       PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf));
589:       PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);
590:       for (l = 0; l < loops_per_scalar; l++) {
591:         DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l);
592:       }
593:       PetscSectionDestroy(&globalSection);
594:       if (newSection) {PetscSectionDestroy(&newSection);}
595:     }
596:   }
597:   /* Cell Fields */
598:   PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL);
599:   if (hasCell || writePartition) {
600:     PetscFPrintf(comm, fp, "CELL_DATA %D\n", totCells);
601:     for (link = vtk->link; link; link = link->next) {
602:       Vec          X = (Vec) link->vec;
603:       PetscSection section = NULL, globalSection;
604:       const char   *name = "";
605:       char         namebuf[256];
606:       PetscInt     enforceDof = PETSC_DETERMINE;

608:       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
609:       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
610:       PetscObjectGetName(link->vec, &name);
611:       PetscObjectQuery(link->vec, "section", (PetscObject*) &section);
612:       if (!section) {
613:         DM           dmX;

615:         VecGetDM(X, &dmX);
616:         if (dmX) {
617:           DMGetLocalSection(dmX, &section);
618:         }
619:       }
620:       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name);
621:       if (link->field >= 0) {
622:         const char *fieldname;

624:         PetscSectionGetFieldName(section, link->field, &fieldname);
625:         PetscSectionGetField(section, link->field, &section);
626:         if (fieldname) {
627:           PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname);
628:         } else {
629:           PetscSNPrintf(namebuf, sizeof(namebuf), "%s%D", name, link->field);
630:         }
631:       } else {
632:         PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name);
633:       }
634:       PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf));
635:       PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);
636:       for (l = 0; l < loops_per_scalar; l++) {
637:         DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l);
638:       }
639:       PetscSectionDestroy(&globalSection);
640:     }
641:     if (writePartition) {
642:       PetscFPrintf(comm, fp, "SCALARS partition int 1\n");
643:       PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");
644:       DMPlexVTKWritePartition_ASCII(dm, fp);
645:     }
646:   }
647:   /* Cleanup */
648:   PetscSectionDestroy(&globalCoordSection);
649:   PetscLayoutDestroy(&vLayout);
650:   PetscFClose(comm, fp);
651:   return(0);
652: }

654: /*@C
655:   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer

657:   Collective

659:   Input Arguments:
660: + odm - The DMPlex specifying the mesh, passed as a PetscObject
661: - viewer - viewer of type VTK

663:   Level: developer

665:   Note:
666:   This function is a callback used by the VTK viewer to actually write the file.
667:   The reason for this odd model is that the VTK file format does not provide any way to write one field at a time.
668:   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.

670: .seealso: PETSCVIEWERVTK
671: @*/
672: PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
673: {
674:   DM             dm = (DM) odm;
675:   PetscBool      isvtk;

681:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);
682:   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
683:   switch (viewer->format) {
684:   case PETSC_VIEWER_ASCII_VTK:
685:     DMPlexVTKWriteAll_ASCII(dm, viewer);
686:     break;
687:   case PETSC_VIEWER_VTK_VTU:
688:     DMPlexVTKWriteAll_VTU(dm, viewer);
689:     break;
690:   default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
691:   }
692:   return(0);
693: }