vtkUnstructuredToFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9  Copyright (C) 2021-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  vtkUnstructuredToFoam
29 
30 Group
31  grpMeshConversionUtilities
32 
33 Description
34  Convert legacy VTK file (ascii) containing an unstructured grid
35  to an OpenFOAM mesh without boundary information.
36 
37 Usage
38  \b vtkUnstructuredToFoam <XXX.vtk>
39 
40  Options:
41  - \par -no-fields
42  Do not attempt to recreate volFields
43 
44 Note
45  The .vtk format does not contain any boundary information.
46  It is purely a description of the internal mesh. This also limits the
47  usefulness of reconstructing the volFields.
48 
49  Not extensively tested.
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #include "argList.H"
54 #include "Time.H"
55 #include "fvMesh.H"
56 #include "IFstream.H"
57 #include "vtkUnstructuredReader.H"
58 
59 #include "columnFvMesh.H"
60 #include "scalarIOField.H"
61 #include "vectorIOField.H"
62 #include "volFields.H"
63 
64 using namespace Foam;
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 template<class Type>
69 void constructVolFields(fvMesh& mesh, const vtkUnstructuredReader& reader)
70 {
71  const auto fields(reader.cellData().csorted<IOField<Type>>());
72  for (const auto& field : fields)
73  {
74  Info<< "Constructing volField " << field.name() << endl;
75 
76  // field is
77  // - cell data followed by
78  // - boundary face data
79 
80 
82  (
83  field.name(),
84  mesh,
85  dimless
86  );
87  auto& fld = tfld.ref();
88  fld.instance() = mesh.time().timeName();
89  fld.writeOpt(IOobject::AUTO_WRITE);
90 
91  // Fill cell values
92  fld.internalFieldRef().field() =
94 
95  // Fill boundary values
96  const auto& map = reader.faceMap();
97  if (map.size())
98  {
99  for (auto& pfld : fld.boundaryFieldRef())
100  {
101  const auto& pp = pfld.patch();
102 
103  forAll(pfld, i)
104  {
105  const label bFacei = pp.patch().offset()+i;
106  pfld[i] = field[map[bFacei]];
107  }
108  }
109  }
110 
111  regIOobject::store(std::move(tfld));
112  }
113 }
114 
115 
116 int main(int argc, char *argv[])
117 {
119  (
120  "Convert legacy VTK file (ascii) containing an unstructured grid"
121  " to an OpenFOAM mesh without boundary information"
122  );
123 
125  argList::addOptionCompat("no-fields", {"noFields", 2106});
126  argList::addArgument("vtk-file", "The input legacy ascii vtk file");
127 
128  #include "setRootCase.H"
129  #include "createTime.H"
130 
131  const bool doFields = !args.found("no-fields");
132 
133  IFstream mshStream(args.get<fileName>(1));
134 
135  vtkUnstructuredReader reader(runTime, mshStream);
136 
137  polyMesh mesh
138  (
139  IOobject
140  (
142  runTime.constant(),
143  runTime
144  ),
145  std::move(reader.points()),
146  reader.cells(),
147  faceListList(),
148  wordList(),
149  wordList(),
150  "defaultFaces",
151  polyPatch::typeName,
152  wordList()
153  );
154 
155  // More precision (for points data)
157 
158  Info<< "Writing mesh ..." << endl;
159 
160  mesh.removeFiles();
161  mesh.write();
162 
163 
164  if (doFields)
165  {
166  // Re-read mesh as fvMesh so we can have fields
167  Info<< "Re-reading mesh ..." << endl;
168  #include "createMesh.H"
169 
170  constructVolFields<scalar>(mesh, reader);
171  constructVolFields<vector>(mesh, reader);
172  constructVolFields<sphericalTensor>(mesh, reader);
173  constructVolFields<symmTensor>(mesh, reader);
174  constructVolFields<tensor>(mesh, reader);
175 
176  // No need to write the mesh, only fields
177  mesh.thisDb().write();
178  }
179 
180 
181  Info<< "End\n" << endl;
182 
183  return 0;
184 }
185 
186 
187 // ************************************************************************* //
List< faceList > faceListList
List of faceList.
Definition: faceListFwd.H:44
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:477
rDeltaTY field()
A class for handling file names.
Definition: fileName.H:72
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
UPtrList< const Type > csorted() const
Return sorted list of objects with a class satisfying isA<Type> or isType<Type> (with Strict) ...
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:519
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
static void noParallel()
Remove the parallel options.
Definition: argList.C:599
Reader for vtk UNSTRUCTURED_GRID legacy files. Supports single CELLS, POINTS etc. entry only...
const pointField & points() const noexcept
Points.
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
static void addOptionCompat(const word &optName, std::pair< const char *, int > compat)
Specify an alias for the option name.
Definition: argList.C:433
const labelList & faceMap() const noexcept
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:299
const cellShapeList & cells() const noexcept
3D cells
const labelList & cellMap() const noexcept
dynamicFvMesh & mesh
void removeFiles(const fileName &instanceDir) const
Remove all files from mesh instance.
Definition: polyMesh.C:1314
const objectRegistry & cellData() const noexcept
Cell based fields.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:406
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
label size() const noexcept
The number of elements in the list.
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition: IOstream.H:459
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1068
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=PatchField< Type >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:131
Input from file stream as an ISstream, normally using std::ifstream for the actual input...
Definition: IFstream.H:51
Required Classes.
List< word > wordList
List of word.
Definition: fileName.H:59
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:52
Automatically write from objectRegistry::writeObject()
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:366
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:188
A primitive field of type <T> with automated input and output.
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.