152 #include "readFields.H" 167 template<
class GeoField>
168 void print(
const char* msg, Ostream& os,
const PtrList<const GeoField>& flds)
175 os <<
' ' << flds[i].name();
182 void print(Ostream& os,
const wordList& flds)
186 os <<
' ' << flds[i];
194 const polyBoundaryMesh& patches,
195 const List<wordRe>& excludePatches
198 DynamicList<label> patchIDs(patches.size());
200 Info<<
"Combining patches:" <<
endl;
204 const polyPatch& pp = patches[
patchi];
208 isA<emptyPolyPatch>(pp)
209 || isA<nonConformalPolyPatch>(pp)
210 || (Pstream::parRun() && isType<processorPolyPatch>(pp))
213 Info<<
" discarding empty/nonConformal/processor patch " 218 Info<<
" excluding patch " <<
patchi 219 <<
" " << pp.name() <<
endl;
224 Info<<
" patch " <<
patchi <<
" " << pp.name() <<
endl;
228 return patchIDs.shrink();
234 int main(
int argc,
char *argv[])
238 "legacy VTK file format writer" 240 timeSelector::addOptions();
247 "only convert the specified fields - eg '(p T U)'" 253 "convert a mesh subset corresponding to the specified cellSet" 259 "restrict conversion to the specified faceSet" 265 "restrict conversion to the specified pointSet" 267 argList::addBoolOption
270 "write in ASCII format instead of binary" 272 argList::addBoolOption
275 "write polyhedral cells without tet/pyramid decomposition" 277 argList::addBoolOption
280 "write surfaceScalarFields (e.g., phi)" 282 argList::addBoolOption
285 "use cell value on patches instead of patch value itself" 287 argList::addBoolOption
290 "do not generate file for mesh, only for patches" 292 argList::addBoolOption
297 argList::addBoolOption
300 "combine all patches into a single file" 306 "a list of patches to exclude - eg '( inlet \".*Wall\" )' " 308 argList::addBoolOption
313 argList::addBoolOption
316 "don't link processor VTK files - parallel only" 318 argList::addBoolOption
321 "use the time name instead of the time index when naming the files" 339 <<
"Using ASCII rather than binary VTK format because " 340 "floatScalar and/or label are not 4 bytes in size." 350 <<
"Using neighbouring cell value instead of patch value" 359 <<
"Outputting cell values only" <<
nl <<
endl;
364 List<wordRe> excludePatches;
369 Info<<
"Not including patches " << excludePatches <<
nl <<
endl;
375 string vtkName = runTime.caseName();
379 vtkName = cellSetName;
381 else if (Pstream::parRun())
386 if (i != string::npos)
388 vtkName = vtkName.substr(i);
401 fileName fvPath(runTime.path()/
"VTK");
404 fileName regionPrefix =
"";
418 || cellSetName.size()
419 || faceSetName.size()
420 || pointSetName.size()
424 Info<<
"Keeping old VTK files in " << fvPath <<
nl <<
endl;
428 Info<<
"Deleting old VTK files in " << fvPath <<
nl <<
endl;
438 vtkMesh vMesh(mesh, cellSetName);
442 HashSet<fileName> allCloudDirs;
445 runTime.setTime(timeDirs[timeI], timeI);
450 runTime.timePath()/regionPrefix/cloud::prefix,
456 IOobjectList sprayObjs
460 cloud::prefix/cloudDirs[i]
463 IOobject* positionsPtr = sprayObjs.lookup(word(
"positions"));
467 if (allCloudDirs.insert(cloudDirs[i]))
469 Info<<
"At time: " << runTime.timeName()
470 <<
" detected cloud directory : " << cloudDirs[i]
480 runTime.setTime(timeDirs[timeI], timeI);
482 Info<<
"Time: " << runTime.timeName() <<
endl;
485 useTimeName ? runTime.timeName() :
Foam::name(runTime.timeIndex());
489 polyMesh::readUpdateState meshState = vMesh.readUpdate();
491 const fvMesh& mesh = vMesh.mesh();
495 meshState == polyMesh::TOPO_CHANGE
496 || meshState == polyMesh::TOPO_PATCH_CHANGE
499 Info<<
" Read new mesh" <<
nl <<
endl;
503 if (faceSetName.size())
506 faceSet
set(
mesh, faceSetName);
511 fileName patchFileName
519 Info<<
" FaceSet : " << patchFileName <<
endl;
527 if (pointSetName.size())
530 pointSet
set(
mesh, pointSetName);
535 fileName patchFileName
543 Info<<
" pointSet : " << patchFileName <<
endl;
552 IOobjectList
objects(mesh, runTime.timeName());
554 HashSet<word> selectedFields;
564 PtrList<const volScalarField::Internal> visf;
565 PtrList<const volVectorField::Internal> vivf;
566 PtrList<const volSphericalTensorField::Internal> visptf;
567 PtrList<const volSymmTensorField::Internal> visytf;
568 PtrList<const volTensorField::Internal> vitf;
570 if (!specifiedFields || selectedFields.size())
573 print(
" volScalarField::Internal :", Info, visf);
576 print(
" volVectorField::Internal :", Info, vivf);
586 print(
" volSphericalTensorFields::Internal :", Info, visptf);
596 print(
" volSymmTensorFields::Internal :", Info, visytf);
599 print(
" volTensorFields::Internal :", Info, vitf);
602 label nVolInternalFields =
612 PtrList<const volScalarField> vsf;
613 PtrList<const volVectorField> vvf;
614 PtrList<const volSphericalTensorField> vsptf;
615 PtrList<const volSymmTensorField> vsytf;
616 PtrList<const volTensorField> vtf;
618 if (!specifiedFields || selectedFields.size())
621 print(
" volScalarFields :", Info, vsf);
624 print(
" volVectorFields :", Info, vvf);
634 print(
" volSphericalTensorFields :", Info, vsptf);
644 print(
" volSymmTensorFields :", Info, vsytf);
647 print(
" volTensorFields :", Info, vtf);
661 Info<<
" pointScalarFields : switched off" 662 <<
" (\"-noPointValues\" (at your option)\n";
663 Info<<
" pointVectorFields : switched off" 664 <<
" (\"-noPointValues\" (at your option)\n";
667 PtrList<const pointScalarField> psf;
668 PtrList<const pointVectorField> pvf;
669 PtrList<const pointSphericalTensorField> psptf;
670 PtrList<const pointSymmTensorField> psytf;
671 PtrList<const pointTensorField> ptf;
673 if (!noPointValues && !(specifiedFields && selectedFields.empty()))
683 print(
" pointScalarFields :", Info, psf);
693 print(
" pointVectorFields :", Info, pvf);
703 print(
" pointSphericalTensorFields :", Info, psptf);
713 print(
" pointSymmTensorFields :", Info, psytf);
723 print(
" pointTensorFields :", Info, ptf);
745 Info<<
" Internal : " << vtkFileName <<
endl;
748 internalWriter writer(vMesh, binary, vtkFileName);
755 1 + nVolInternalFields + nVolFields
759 writer.writeCellIDs();
764 writer.write(visptf);
765 writer.write(visytf);
780 vMesh.nFieldPoints(),
781 nVolFields + nPointFields
792 volPointInterpolation pInterp(mesh);
793 writer.write(pInterp, vsf);
794 writer.write(pInterp, vvf);
795 writer.write(pInterp, vsptf);
796 writer.write(pInterp, vsytf);
797 writer.write(pInterp, vtf);
809 PtrList<const surfaceScalarField> ssf;
818 print(
" surfScalarFields :", Info, ssf);
820 PtrList<const surfaceVectorField> svf;
829 print(
" surfVectorFields :", Info, svf);
831 if (ssf.size() + svf.size() > 0)
834 label sz = svf.size();
836 svf.setSize(sz + ssf.size());
844 svf.set(sz+i, ssfiPtr);
848 mkDir(fvPath /
"surfaceFields");
850 fileName surfFileName
877 const polyBoundaryMesh& patches = mesh.boundaryMesh();
879 const labelList patchIDs(getSelectedPatches(patches, excludePatches));
883 mkDir(fvPath/
"allPatches");
885 fileName patchFileName;
887 if (vMesh.useSubMesh())
890 fvPath/
"allPatches"/cellSetName
898 fvPath/
"allPatches"/
"allPatches" 904 Info<<
" Combined patches : " << patchFileName <<
endl;
912 getSelectedPatches(patches, excludePatches)
924 writer.writePatchIDs();
954 const polyPatch& pp = patches[patchIDs[i]];
956 mkDir(fvPath/pp.name());
958 fileName patchFileName;
960 if (vMesh.useSubMesh())
963 fvPath/pp.name()/cellSetName
971 fvPath/pp.name()/pp.name()
977 Info<<
" Patch : " << patchFileName <<
endl;
997 writer.writePatchIDs();
1002 writer.write(vsptf);
1003 writer.write(vsytf);
1019 writer.write(psptf);
1020 writer.write(psytf);
1023 PrimitivePatchInterpolation<primitivePatch> pInter
1029 writer.write(pInter, vsf);
1030 writer.write(pInter, vvf);
1031 writer.write(pInter, vsptf);
1032 writer.write(pInter, vsytf);
1033 writer.write(pInter, vtf);
1046 PtrList<const surfaceScalarField> ssf;
1055 print(
" surfScalarFields :", Info, ssf);
1057 PtrList<const surfaceVectorField> svf;
1066 print(
" surfVectorFields :", Info, svf);
1072 const faceZone& fz = zones[zoneI];
1074 mkDir(fvPath/fz.name());
1076 fileName patchFileName;
1078 if (vMesh.useSubMesh())
1081 fvPath/fz.name()/cellSetName
1089 fvPath/fz.name()/fz.name()
1095 Info<<
" FaceZone : " << patchFileName <<
endl;
1099 IndirectList<face>(mesh.faces(), fz),
1103 surfaceMeshWriter writer
1116 ssf.size() + svf.size()
1134 const fileName& cloudName = iter.key();
1137 mkDir(fvPath/cloud::prefix/cloudName);
1139 fileName lagrFileName
1141 fvPath/cloud::prefix/cloudName/cloudName
1142 +
"_" + timeDesc +
".vtk" 1145 Info<<
" Lagrangian: " << lagrFileName <<
endl;
1148 IOobjectList sprayObjs
1155 IOobject* positionsPtr = sprayObjs.lookup(word(
"positions"));
1159 wordList labelNames(sprayObjs.names(labelIOField::typeName));
1161 print(Info, labelNames);
1163 wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1164 Info<<
" scalars :";
1165 print(Info, scalarNames);
1167 wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1168 Info<<
" vectors :";
1169 print(Info, vectorNames);
1175 sphericalTensorIOField::typeName
1178 Info<<
" spherical tensors :";
1179 print(Info, sphereNames);
1185 symmTensorIOField::typeName
1188 Info<<
" symm tensors :";
1189 print(Info, symmNames);
1191 wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
1192 Info<<
" tensors :";
1193 print(Info, tensorNames);
1195 lagrangianWriter writer
1205 writer.writeFieldsHeader
1208 + scalarNames.size()
1209 + vectorNames.size()
1210 + sphereNames.size()
1212 + tensorNames.size()
1216 writer.writeIOField<
label>(labelNames);
1217 writer.writeIOField<scalar>(scalarNames);
1218 writer.writeIOField<
vector>(vectorNames);
1221 writer.writeIOField<
tensor>(tensorNames);
1225 lagrangianWriter writer
1235 writer.writeFieldsHeader(0);
1247 if (Pstream::parRun() && doLinks)
1249 mkDir(runTime.globalPath()/
"VTK");
1250 chDir(runTime.globalPath()/
"VTK");
1252 Info<<
"Linking all processor files to " << runTime.globalPath()/
"VTK" 1259 /
"processor" +
name(Pstream::myProcNo())
1264 label sz = dirs.size();
1265 dirs.setSize(sz + 1);
1274 fileName procFile(procVTK/dirs[i]/subFiles[j]);
1282 +
name(Pstream::myProcNo())
1291 Info<<
"End\n" <<
endl;
List< instant > instantList
List of instants.
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
#define forAll(list, i)
Loop across all elements in list.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
bool optionFound(const word &opt) const
Return true if the named option is found.
void writePointSet(const bool binary, const vtkMesh &vMesh, const pointSet &set, const fileName &fileName)
bool chDir(const fileName &dir)
Change the current directory to the one given and return true,.
void writeCellDataHeader(std::ostream &, const label nCells, const label nFields)
Vector< scalar > vector
A scalar version of the templated Vector.
Operations on lists of strings.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
void writeSurfFields(const bool binary, const vtkMesh &vMesh, const fileName &fileName, const UPtrList< const surfaceVectorField > &surfVectorFields)
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
float floatScalar
Float precision floating point scalar type.
Write faceSet to vtk polydata file. Only one data which is original faceID.
void writeFaceSet(const bool binary, const vtkMesh &vMesh, const faceSet &set, const fileName &fileName)
List< label > labelList
A List of labels.
fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true)
Read a directory and return the entries as a string list.
graph_traits< Graph >::vertices_size_type size_type
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
bool rmDir(const fileName &)
Remove a directory and its contents.
Write a patch with its data.
virtual void rename(const word &newName)
Rename.
Write pointSet to vtk polydata file. Only one data which is original pointID.
MeshZones< faceZone, polyMesh > meshFaceZones
A MeshZones with the type faceZone.
word name(const complex &)
Return a string representation of a complex.
List< word > wordList
A List of words.
#define WarningInFunction
Report a warning using Foam::Warning.
void writePointDataHeader(std::ostream &, const label nPoints, const label nFields)
const word cloudName(propsDict.lookup("cloudName"))
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars.
Foam::argList args(argc, argv)
List< fileName > fileNameList
A List of fileNames.
Tensor< scalar > tensor
Tensor of scalars.
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.