foamToTecplot360.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  foamToTecplot360
26 
27 Description
28  Tecplot binary file format writer.
29 
30 Usage
31  \b foamToTecplot360 [OPTION]
32 
33  Options:
34  - \par -fields <names>
35  Convert selected fields only. For example,
36  \verbatim
37  -fields '( p T U )'
38  \endverbatim
39  The quoting is required to avoid shell expansions and to pass the
40  information as a single argument.
41 
42  - \par -cellSet <name>
43 
44  - \par -faceSet <name>
45  Restrict conversion to the cellSet, faceSet.
46 
47  - \par -nearCellValue
48  Output cell value on patches instead of patch value itself
49 
50  - \par -noInternal
51  Do not generate file for mesh, only for patches
52 
53  - \par -noPointValues
54  No pointFields
55 
56  - \par -noFaceZones
57  No faceZones
58 
59  - \par -excludePatches <patchNames>
60  Specify patches (wildcards) to exclude. For example,
61  \verbatim
62  -excludePatches '( inlet_1 inlet_2 "proc.*")'
63  \endverbatim
64  The quoting is required to avoid shell expansions and to pass the
65  information as a single argument. The double quotes denote a regular
66  expression.
67 
68  - \par -useTimeName
69  use the time index in the VTK file name instead of the time index
70 
71 \*---------------------------------------------------------------------------*/
72 
73 #include "argList.H"
74 #include "pointMesh.H"
75 #include "volPointInterpolation.H"
76 #include "emptyPolyPatch.H"
77 #include "labelIOField.H"
78 #include "scalarIOField.H"
79 #include "sphericalTensorIOField.H"
80 #include "symmTensorIOField.H"
81 #include "tensorIOField.H"
82 #include "passiveParticleCloud.H"
83 #include "faceSet.H"
84 #include "stringListOps.H"
85 #include "wordRe.H"
86 
87 #include "vtkMesh.H"
88 #include "readFields.H"
89 #include "tecplotWriter.H"
90 
91 #include "TECIO.h"
92 
93 // Note: needs to be after TECIO to prevent Foam::Time conflicting with
94 // Xlib Time.
95 #include "Time.H"
96 
97 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
98 
99 template<class GeoField>
100 void print(const char* msg, Ostream& os, const PtrList<GeoField>& flds)
101 {
102  if (flds.size())
103  {
104  os << msg;
105  forAll(flds, i)
106  {
107  os << ' ' << flds[i].name();
108  }
109  os << endl;
110  }
111 }
112 
113 
114 void print(Ostream& os, const wordList& flds)
115 {
116  forAll(flds, i)
117  {
118  os << ' ' << flds[i];
119  }
120  os << endl;
121 }
122 
123 
124 labelList getSelectedPatches
125 (
126  const polyBoundaryMesh& patches,
127  const List<wordRe>& excludePatches // HashSet<word>& excludePatches
128 )
129 {
130  DynamicList<label> patchIDs(patches.size());
131 
132  Info<< "Combining patches:" << endl;
133 
135  {
136  const polyPatch& pp = patches[patchi];
137 
138  if
139  (
140  isType<emptyPolyPatch>(pp)
141  || (Pstream::parRun() && isType<processorPolyPatch>(pp))
142  )
143  {
144  Info<< " discarding empty/processor patch " << patchi
145  << " " << pp.name() << endl;
146  }
147  else if (findStrings(excludePatches, pp.name()))
148  {
149  Info<< " excluding patch " << patchi
150  << " " << pp.name() << endl;
151  }
152  else
153  {
154  patchIDs.append(patchi);
155  Info<< " patch " << patchi << " " << pp.name() << endl;
156  }
157  }
158  return patchIDs.shrink();
159 }
160 
161 
162 
163 
164 
165 int main(int argc, char *argv[])
166 {
167  argList::addNote
168  (
169  "Tecplot binary file format writer"
170  );
171 
172  timeSelector::addOptions();
173  #include "addRegionOption.H"
174 
175  argList::addOption
176  (
177  "fields",
178  "names",
179  "convert selected fields only. eg, '(p T U)'"
180  );
181  argList::addOption
182  (
183  "cellSet",
184  "name",
185  "restrict conversion to the specified cellSet"
186  );
187  argList::addOption
188  (
189  "faceSet",
190  "name",
191  "restrict conversion to the specified cellSet"
192  );
193  argList::addBoolOption
194  (
195  "nearCellValue",
196  "output cell value on patches instead of patch value itself"
197  );
198  argList::addBoolOption
199  (
200  "noInternal",
201  "do not generate file for mesh, only for patches"
202  );
203  argList::addBoolOption
204  (
205  "noPointValues",
206  "no pointFields"
207  );
208  argList::addOption
209  (
210  "excludePatches",
211  "patches (wildcards) to exclude"
212  );
213  argList::addBoolOption
214  (
215  "noFaceZones",
216  "no faceZones"
217  );
218 
219  #include "setRootCase.H"
220  #include "createTime.H"
221 
222  const bool doWriteInternal = !args.optionFound("noInternal");
223  const bool doFaceZones = !args.optionFound("noFaceZones");
224  const bool nearCellValue = args.optionFound("nearCellValue");
225  const bool noPointValues = args.optionFound("noPointValues");
226 
227  if (nearCellValue)
228  {
230  << "Using neighbouring cell value instead of patch value"
231  << nl << endl;
232  }
233 
234  if (noPointValues)
235  {
237  << "Outputting cell values only" << nl << endl;
238  }
239 
240  List<wordRe> excludePatches;
241  if (args.optionFound("excludePatches"))
242  {
243  args.optionLookup("excludePatches")() >> excludePatches;
244 
245  Info<< "Not including patches " << excludePatches << nl << endl;
246  }
247 
248  word cellSetName;
249  string vtkName;
250 
251  if (args.optionReadIfPresent("cellSet", cellSetName))
252  {
253  vtkName = cellSetName;
254  }
255  else if (Pstream::parRun())
256  {
257  // Strip off leading casename, leaving just processor_DDD ending.
258  vtkName = runTime.caseName();
259 
260  string::size_type i = vtkName.rfind("processor");
261 
262  if (i != string::npos)
263  {
264  vtkName = vtkName.substr(i);
265  }
266  }
267  else
268  {
269  vtkName = runTime.caseName();
270  }
271 
272 
273  const instantList timeDirs = timeSelector::select0(runTime, args);
274 
276 
277  // TecplotData/ directory in the case
278  fileName fvPath(runTime.path()/"Tecplot360");
279  // Directory of mesh (region0 gets filtered out)
280  fileName regionPrefix = "";
281 
282  if (regionName != polyMesh::defaultRegion)
283  {
284  fvPath = fvPath/regionName;
285  regionPrefix = regionName;
286  }
287 
288  if (isDir(fvPath))
289  {
290  if
291  (
292  args.optionFound("time")
293  || args.optionFound("latestTime")
294  || cellSetName.size()
295  || regionName != polyMesh::defaultRegion
296  )
297  {
298  Info<< "Keeping old files in " << fvPath << nl << endl;
299  }
300  else
301  {
302  Info<< "Deleting old VTK files in " << fvPath << nl << endl;
303 
304  rmDir(fvPath);
305  }
306  }
307 
308  mkDir(fvPath);
309 
310 
311  // mesh wrapper; does subsetting and decomposition
312  vtkMesh vMesh(mesh, cellSetName);
313 
314  forAll(timeDirs, timeI)
315  {
316  runTime.setTime(timeDirs[timeI], timeI);
317 
318  Info<< "Time: " << runTime.name() << endl;
319 
320  const word timeDesc = name(timeI); // name(runTime.timeIndex());
321 
322  // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
323  // decomposition.
324  fvMesh::readUpdateState meshState = vMesh.readUpdate();
325 
326  const fvMesh& mesh = vMesh.mesh();
327 
328  INTEGER4 nFaceNodes = 0;
329  forAll(mesh.faces(), facei)
330  {
331  nFaceNodes += mesh.faces()[facei].size();
332  }
333 
334 
335  // Read all fields on the new mesh
336  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
337 
338  // Search for list of objects for this time
339  IOobjectList objects(mesh, runTime.name());
340 
341  HashSet<word> selectedFields;
342  if (args.optionFound("fields"))
343  {
344  args.optionLookup("fields")() >> selectedFields;
345  }
346 
347  // Construct the vol fields (on the original mesh if subsetted)
348 
349  PtrList<volScalarField> vsf;
350  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vsf);
351  print(" volScalarFields :", Info, vsf);
352 
353  PtrList<volVectorField> vvf;
354  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vvf);
355  print(" volVectorFields :", Info, vvf);
356 
357  PtrList<volSphericalTensorField> vSpheretf;
358  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSpheretf);
359  print(" volSphericalTensorFields :", Info, vSpheretf);
360 
361  PtrList<volSymmTensorField> vSymmtf;
362  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSymmtf);
363  print(" volSymmTensorFields :", Info, vSymmtf);
364 
365  PtrList<volTensorField> vtf;
366  readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vtf);
367  print(" volTensorFields :", Info, vtf);
368 
369 
370  // Construct pointMesh only if necessary since constructs edge
371  // addressing (expensive on polyhedral meshes)
372  if (noPointValues)
373  {
374  Info<< " pointScalarFields : switched off"
375  << " (\"-noPointValues\" (at your option)\n";
376  Info<< " pointVectorFields : switched off"
377  << " (\"-noPointValues\" (at your option)\n";
378  }
379 
380  PtrList<pointScalarField> psf;
381  PtrList<pointVectorField> pvf;
382  // PtrList<pointSphericalTensorField> pSpheretf;
383  // PtrList<pointSymmTensorField> pSymmtf;
384  // PtrList<pointTensorField> ptf;
385 
386 
387  if (!noPointValues)
388  {
390  // const volPointInterpolation& pInterp = volPointInterpolation::New
391  //(
392  // mesh
393  //);
394  //
395  // label nPsf = psf.size();
396  // psf.setSize(nPsf+vsf.size());
397  // forAll(vsf, i)
398  //{
399  // Info<< "Interpolating " << vsf[i].name() << endl;
400  // tmp<pointScalarField> tvsf(pInterp.interpolate(vsf[i]));
401  // tvsf().rename(vsf[i].name() + "_point");
402  // psf.set(nPsf, tvsf);
403  // nPsf++;
404  //}
405  //
406  // label nPvf = pvf.size();
407  // pvf.setSize(vvf.size());
408  // forAll(vvf, i)
409  //{
410  // Info<< "Interpolating " << vvf[i].name() << endl;
411  // tmp<pointVectorField> tvvf(pInterp.interpolate(vvf[i]));
412  // tvvf().rename(vvf[i].name() + "_point");
413  // pvf.set(nPvf, tvvf);
414  // nPvf++;
415  //}
416 
417  readFields
418  (
419  vMesh,
420  pointMesh::New(vMesh.baseMesh()),
421  objects,
422  selectedFields,
423  psf
424  );
425  print(" pointScalarFields :", Info, psf);
426 
427  readFields
428  (
429  vMesh,
430  pointMesh::New(vMesh.baseMesh()),
431  objects,
432  selectedFields,
433  pvf
434  );
435  print(" pointVectorFields :", Info, pvf);
436 
437  // readFields
438  //(
439  // vMesh,
440  // pointMesh::New(vMesh.baseMesh()),
441  // objects,
442  // selectedFields,
443  // pSpheretf
444  //);
445  // print(" pointSphericalTensorFields :", Info, pSpheretf);
446  //
447  // readFields
448  //(
449  // vMesh,
450  // pointMesh::New(vMesh.baseMesh()),
451  // objects,
452  // selectedFields,
453  // pSymmtf
454  //);
455  // print(" pointSymmTensorFields :", Info, pSymmtf);
456  //
457  // readFields
458  //(
459  // vMesh,
460  // pointMesh::New(vMesh.baseMesh()),
461  // objects,
462  // selectedFields,
463  // ptf
464  //);
465  // print(" pointTensorFields :", Info, ptf);
466  }
467  Info<< endl;
468 
469 
470  // Get field names
471  // ~~~~~~~~~~~~~~~
472 
473  string varNames;
474  DynamicList<INTEGER4> varLocation;
475 
476  string cellVarNames;
477  DynamicList<INTEGER4> cellVarLocation;
478 
479  // volFields
480  tecplotWriter::getTecplotNames
481  (
482  vsf,
483  ValueLocation_CellCentered,
484  varNames,
485  varLocation
486  );
487  tecplotWriter::getTecplotNames
488  (
489  vsf,
490  ValueLocation_CellCentered,
491  cellVarNames,
492  cellVarLocation
493  );
494 
495  tecplotWriter::getTecplotNames
496  (
497  vvf,
498  ValueLocation_CellCentered,
499  varNames,
500  varLocation
501  );
502  tecplotWriter::getTecplotNames
503  (
504  vvf,
505  ValueLocation_CellCentered,
506  cellVarNames,
507  cellVarLocation
508  );
509 
510  tecplotWriter::getTecplotNames
511  (
512  vSpheretf,
513  ValueLocation_CellCentered,
514  varNames,
515  varLocation
516  );
517  tecplotWriter::getTecplotNames
518  (
519  vSpheretf,
520  ValueLocation_CellCentered,
521  cellVarNames,
522  cellVarLocation
523  );
524 
525  tecplotWriter::getTecplotNames
526  (
527  vSymmtf,
528  ValueLocation_CellCentered,
529  varNames,
530  varLocation
531  );
532  tecplotWriter::getTecplotNames
533  (
534  vSymmtf,
535  ValueLocation_CellCentered,
536  cellVarNames,
537  cellVarLocation
538  );
539 
540  tecplotWriter::getTecplotNames
541  (
542  vtf,
543  ValueLocation_CellCentered,
544  varNames,
545  varLocation
546  );
547  tecplotWriter::getTecplotNames
548  (
549  vtf,
550  ValueLocation_CellCentered,
551  cellVarNames,
552  cellVarLocation
553  );
554 
555 
556 
557  // pointFields
558  tecplotWriter::getTecplotNames
559  (
560  psf,
561  ValueLocation_Nodal,
562  varNames,
563  varLocation
564  );
565 
566  tecplotWriter::getTecplotNames
567  (
568  pvf,
569  ValueLocation_Nodal,
570  varNames,
571  varLocation
572  );
573 
574  // strandID (= piece id. Gets incremented for every piece of geometry
575  // that is output)
576  INTEGER4 strandID = 1;
577 
578 
579  if (meshState != fvMesh::UNCHANGED)
580  {
581  if (doWriteInternal)
582  {
583  // Output mesh and fields
584  fileName vtkFileName
585  (
586  fvPath/vtkName
587  + "_"
588  + timeDesc
589  + ".plt"
590  );
591 
592  tecplotWriter writer(runTime);
593 
594  string allVarNames = string("X Y Z ") + varNames;
595  DynamicList<INTEGER4> allVarLocation;
596  allVarLocation.append(ValueLocation_Nodal);
597  allVarLocation.append(ValueLocation_Nodal);
598  allVarLocation.append(ValueLocation_Nodal);
599  allVarLocation.append(varLocation);
600 
601  writer.writeInit
602  (
603  runTime.caseName(),
604  allVarNames,
605  vtkFileName,
606  DataFileType_Full
607  );
608 
609  writer.writePolyhedralZone
610  (
611  mesh.name(), // regionName
612  strandID++, // strandID
613  mesh,
614  allVarLocation,
615  nFaceNodes
616  );
617 
618  // Write coordinates
619  writer.writeField(mesh.points().component(0)());
620  writer.writeField(mesh.points().component(1)());
621  writer.writeField(mesh.points().component(2)());
622 
623  // Write all fields
624  forAll(vsf, i)
625  {
626  writer.writeField(vsf[i]);
627  }
628  forAll(vvf, i)
629  {
630  writer.writeField(vvf[i]);
631  }
632  forAll(vSpheretf, i)
633  {
634  writer.writeField(vSpheretf[i]);
635  }
636  forAll(vSymmtf, i)
637  {
638  writer.writeField(vSymmtf[i]);
639  }
640  forAll(vtf, i)
641  {
642  writer.writeField(vtf[i]);
643  }
644 
645  forAll(psf, i)
646  {
647  writer.writeField(psf[i]);
648  }
649  forAll(pvf, i)
650  {
651  writer.writeField(pvf[i]);
652  }
653 
654  writer.writeConnectivity(mesh);
655  writer.writeEnd();
656  }
657  }
658  else
659  {
660  if (doWriteInternal)
661  {
662  if (timeI == 0)
663  {
664  // Output static mesh only
665  fileName vtkFileName
666  (
667  fvPath/vtkName
668  + "_grid_"
669  + timeDesc
670  + ".plt"
671  );
672 
673  tecplotWriter writer(runTime);
674 
675  writer.writeInit
676  (
677  runTime.caseName(),
678  "X Y Z",
679  vtkFileName,
680  DataFileType_Grid
681  );
682 
683  writer.writePolyhedralZone
684  (
685  mesh.name(), // regionName
686  strandID, // strandID
687  mesh,
688  List<INTEGER4>(3, ValueLocation_Nodal),
689  nFaceNodes
690  );
691 
692  // Write coordinates
693  writer.writeField(mesh.points().component(0)());
694  writer.writeField(mesh.points().component(1)());
695  writer.writeField(mesh.points().component(2)());
696 
697  writer.writeConnectivity(mesh);
698  writer.writeEnd();
699  }
700 
701  // Output solution file
702  fileName vtkFileName
703  (
704  fvPath/vtkName
705  + "_"
706  + timeDesc
707  + ".plt"
708  );
709 
710  tecplotWriter writer(runTime);
711 
712  writer.writeInit
713  (
714  runTime.caseName(),
715  varNames,
716  vtkFileName,
717  DataFileType_Solution
718  );
719 
720  writer.writePolyhedralZone
721  (
722  mesh.name(), // regionName
723  strandID++, // strandID
724  mesh,
725  varLocation,
726  0
727  );
728 
729  // Write all fields
730  forAll(vsf, i)
731  {
732  writer.writeField(vsf[i]);
733  }
734  forAll(vvf, i)
735  {
736  writer.writeField(vvf[i]);
737  }
738  forAll(vSpheretf, i)
739  {
740  writer.writeField(vSpheretf[i]);
741  }
742  forAll(vSymmtf, i)
743  {
744  writer.writeField(vSymmtf[i]);
745  }
746  forAll(vtf, i)
747  {
748  writer.writeField(vtf[i]);
749  }
750 
751  forAll(psf, i)
752  {
753  writer.writeField(psf[i]);
754  }
755  forAll(pvf, i)
756  {
757  writer.writeField(pvf[i]);
758  }
759  writer.writeEnd();
760  }
761  }
762 
763 
764  //---------------------------------------------------------------------
765  //
766  // Write faceSet
767  //
768  //---------------------------------------------------------------------
769 
770  if (args.optionFound("faceSet"))
771  {
772  // Load the faceSet
773  const word setName = args["faceSet"];
774  labelList faceLabels(faceSet(mesh, setName).toc());
775 
776  // Filename as if patch with same name.
777  mkDir(fvPath/setName);
778 
779  fileName patchFileName
780  (
781  fvPath/setName/setName
782  + "_"
783  + timeDesc
784  + ".plt"
785  );
786 
787  Info<< " FaceSet : " << patchFileName << endl;
788 
789  tecplotWriter writer(runTime);
790 
791  string allVarNames = string("X Y Z ") + cellVarNames;
792  DynamicList<INTEGER4> allVarLocation;
793  allVarLocation.append(ValueLocation_Nodal);
794  allVarLocation.append(ValueLocation_Nodal);
795  allVarLocation.append(ValueLocation_Nodal);
796  allVarLocation.append(cellVarLocation);
797 
798  writer.writeInit
799  (
800  runTime.caseName(),
801  cellVarNames,
802  patchFileName,
803  DataFileType_Full
804  );
805 
806  const indirectPrimitivePatch ipp
807  (
808  IndirectList<face>(mesh.faces(), faceLabels),
809  mesh.points()
810  );
811 
812  writer.writePolygonalZone
813  (
814  setName,
815  strandID++,
816  ipp,
817  allVarLocation
818  );
819 
820  // Write coordinates
821  writer.writeField(ipp.localPoints().component(0)());
822  writer.writeField(ipp.localPoints().component(1)());
823  writer.writeField(ipp.localPoints().component(2)());
824 
825  // Write all volfields
826  forAll(vsf, i)
827  {
828  writer.writeField
829  (
830  writer.getFaceField
831  (
832  linearInterpolate(vsf[i])(),
833  faceLabels
834  )()
835  );
836  }
837  forAll(vvf, i)
838  {
839  writer.writeField
840  (
841  writer.getFaceField
842  (
843  linearInterpolate(vvf[i])(),
844  faceLabels
845  )()
846  );
847  }
848  forAll(vSpheretf, i)
849  {
850  writer.writeField
851  (
852  writer.getFaceField
853  (
854  linearInterpolate(vSpheretf[i])(),
855  faceLabels
856  )()
857  );
858  }
859  forAll(vSymmtf, i)
860  {
861  writer.writeField
862  (
863  writer.getFaceField
864  (
865  linearInterpolate(vSymmtf[i])(),
866  faceLabels
867  )()
868  );
869  }
870  forAll(vtf, i)
871  {
872  writer.writeField
873  (
874  writer.getFaceField
875  (
876  linearInterpolate(vtf[i])(),
877  faceLabels
878  )()
879  );
880  }
881  writer.writeConnectivity(ipp);
882 
883  continue;
884  }
885 
886 
887 
888  //---------------------------------------------------------------------
889  //
890  // Write patches as multi-zone file
891  //
892  //---------------------------------------------------------------------
893 
894  const polyBoundaryMesh& patches = mesh.boundaryMesh();
895 
896  labelList patchIDs(getSelectedPatches(patches, excludePatches));
897 
898  mkDir(fvPath/"boundaryMesh");
899 
900  fileName patchFileName;
901 
902  if (vMesh.useSubMesh())
903  {
904  patchFileName =
905  fvPath/"boundaryMesh"/cellSetName
906  + "_"
907  + timeDesc
908  + ".plt";
909  }
910  else
911  {
912  patchFileName =
913  fvPath/"boundaryMesh"/"boundaryMesh"
914  + "_"
915  + timeDesc
916  + ".plt";
917  }
918 
919  Info<< " Combined patches : " << patchFileName << endl;
920 
921  tecplotWriter writer(runTime);
922 
923  string allVarNames = string("X Y Z ") + varNames;
924  DynamicList<INTEGER4> allVarLocation;
925  allVarLocation.append(ValueLocation_Nodal);
926  allVarLocation.append(ValueLocation_Nodal);
927  allVarLocation.append(ValueLocation_Nodal);
928  allVarLocation.append(varLocation);
929 
930  writer.writeInit
931  (
932  runTime.caseName(),
933  allVarNames,
934  patchFileName,
935  DataFileType_Full
936  );
937 
938  forAll(patchIDs, i)
939  {
940  label patchID = patchIDs[i];
941  const polyPatch& pp = patches[patchID];
942  // INTEGER4 strandID = 1 + i;
943 
944  if (pp.size() > 0)
945  {
946  Info<< " Writing patch " << patchID << "\t" << pp.name()
947  << "\tstrand:" << strandID << nl << endl;
948 
949  const indirectPrimitivePatch ipp
950  (
951  IndirectList<face>(pp, identityMap(pp.size())),
952  pp.points()
953  );
954 
955  writer.writePolygonalZone
956  (
957  pp.name(),
958  strandID++, // strandID,
959  ipp,
960  allVarLocation
961  );
962 
963  // Write coordinates
964  writer.writeField(ipp.localPoints().component(0)());
965  writer.writeField(ipp.localPoints().component(1)());
966  writer.writeField(ipp.localPoints().component(2)());
967 
968  // Write all fields
969  forAll(vsf, i)
970  {
971  writer.writeField
972  (
973  writer.getPatchField
974  (
975  nearCellValue,
976  vsf[i],
977  patchID
978  )()
979  );
980  }
981  forAll(vvf, i)
982  {
983  writer.writeField
984  (
985  writer.getPatchField
986  (
987  nearCellValue,
988  vvf[i],
989  patchID
990  )()
991  );
992  }
993  forAll(vSpheretf, i)
994  {
995  writer.writeField
996  (
997  writer.getPatchField
998  (
999  nearCellValue,
1000  vSpheretf[i],
1001  patchID
1002  )()
1003  );
1004  }
1005  forAll(vSymmtf, i)
1006  {
1007  writer.writeField
1008  (
1009  writer.getPatchField
1010  (
1011  nearCellValue,
1012  vSymmtf[i],
1013  patchID
1014  )()
1015  );
1016  }
1017  forAll(vtf, i)
1018  {
1019  writer.writeField
1020  (
1021  writer.getPatchField
1022  (
1023  nearCellValue,
1024  vtf[i],
1025  patchID
1026  )()
1027  );
1028  }
1029 
1030  forAll(psf, i)
1031  {
1032  writer.writeField
1033  (
1034  psf[i].boundaryField()[patchID].patchInternalField()()
1035  );
1036  }
1037  forAll(pvf, i)
1038  {
1039  writer.writeField
1040  (
1041  pvf[i].boundaryField()[patchID].patchInternalField()()
1042  );
1043  }
1044 
1045  writer.writeConnectivity(ipp);
1046  }
1047  else
1048  {
1049  Info<< " Skipping zero sized patch " << patchID
1050  << "\t" << pp.name()
1051  << nl << endl;
1052  }
1053  }
1054  writer.writeEnd();
1055  Info<< endl;
1056 
1057 
1058 
1059  //---------------------------------------------------------------------
1060  //
1061  // Write faceZones as multi-zone file
1062  //
1063  //---------------------------------------------------------------------
1064 
1065  const faceZoneList& zones = mesh.faceZones();
1066 
1067  if (doFaceZones && zones.size() > 0)
1068  {
1069  mkDir(fvPath/"faceZones");
1070 
1071  fileName patchFileName;
1072 
1073  if (vMesh.useSubMesh())
1074  {
1075  patchFileName =
1076  fvPath/"faceZones"/cellSetName
1077  + "_"
1078  + timeDesc
1079  + ".plt";
1080  }
1081  else
1082  {
1083  patchFileName =
1084  fvPath/"faceZones"/"faceZones"
1085  + "_"
1086  + timeDesc
1087  + ".plt";
1088  }
1089 
1090  Info<< " FaceZone : " << patchFileName << endl;
1091 
1092  tecplotWriter writer(runTime);
1093 
1094  string allVarNames = string("X Y Z ") + cellVarNames;
1095  DynamicList<INTEGER4> allVarLocation;
1096  allVarLocation.append(ValueLocation_Nodal);
1097  allVarLocation.append(ValueLocation_Nodal);
1098  allVarLocation.append(ValueLocation_Nodal);
1099  allVarLocation.append(cellVarLocation);
1100 
1101  writer.writeInit
1102  (
1103  runTime.caseName(),
1104  allVarNames,
1105  patchFileName,
1106  DataFileType_Full
1107  );
1108 
1109  forAll(zones, zoneI)
1110  {
1111  const faceZone& pp = zones[zoneI];
1112 
1113  if (pp.size() > 0)
1114  {
1115  const indirectPrimitivePatch ipp
1116  (
1117  IndirectList<face>(mesh.faces(), pp),
1118  mesh.points()
1119  );
1120 
1121  writer.writePolygonalZone
1122  (
1123  pp.name(),
1124  strandID++, //1+patchIDs.size()+zoneI, // strandID,
1125  ipp,
1126  allVarLocation
1127  );
1128 
1129  // Write coordinates
1130  writer.writeField(ipp.localPoints().component(0)());
1131  writer.writeField(ipp.localPoints().component(1)());
1132  writer.writeField(ipp.localPoints().component(2)());
1133 
1134  // Write all volfields
1135  forAll(vsf, i)
1136  {
1137  writer.writeField
1138  (
1139  writer.getFaceField
1140  (
1141  linearInterpolate(vsf[i])(),
1142  pp
1143  )()
1144  );
1145  }
1146  forAll(vvf, i)
1147  {
1148  writer.writeField
1149  (
1150  writer.getFaceField
1151  (
1152  linearInterpolate(vvf[i])(),
1153  pp
1154  )()
1155  );
1156  }
1157  forAll(vSpheretf, i)
1158  {
1159  writer.writeField
1160  (
1161  writer.getFaceField
1162  (
1163  linearInterpolate(vSpheretf[i])(),
1164  pp
1165  )()
1166  );
1167  }
1168  forAll(vSymmtf, i)
1169  {
1170  writer.writeField
1171  (
1172  writer.getFaceField
1173  (
1174  linearInterpolate(vSymmtf[i])(),
1175  pp
1176  )()
1177  );
1178  }
1179  forAll(vtf, i)
1180  {
1181  writer.writeField
1182  (
1183  writer.getFaceField
1184  (
1185  linearInterpolate(vtf[i])(),
1186  pp
1187  )()
1188  );
1189  }
1190 
1191  writer.writeConnectivity(ipp);
1192  }
1193  else
1194  {
1195  Info<< " Skipping zero sized faceZone " << zoneI
1196  << "\t" << pp.name()
1197  << nl << endl;
1198  }
1199  }
1200  writer.writeEnd();
1201  Info<< endl;
1202  }
1203 
1204 
1205 
1206  //---------------------------------------------------------------------
1207  //
1208  // Write lagrangian data
1209  //
1210  //---------------------------------------------------------------------
1211 
1212  fileNameList cloudDirs
1213  (
1214  readDir
1215  (
1216  runTime.timePath()/regionPrefix/lagrangian::cloud::prefix,
1217  fileType::directory
1218  )
1219  );
1220 
1221  forAll(cloudDirs, cloudI)
1222  {
1223  IOobjectList sprayObjs
1224  (
1225  mesh,
1226  runTime.name(),
1227  lagrangian::cloud::prefix/cloudDirs[cloudI]
1228  );
1229 
1230  IOobject* positionsPtr = sprayObjs.lookup("positions");
1231 
1232  if (positionsPtr)
1233  {
1234  mkDir(fvPath/lagrangian::cloud::prefix/cloudDirs[cloudI]);
1235 
1236  fileName lagrFileName
1237  (
1238  fvPath
1239  /lagrangian::cloud::prefix
1240  /cloudDirs[cloudI]
1241  /cloudDirs[cloudI]
1242  + "_" + timeDesc + ".plt"
1243  );
1244 
1245  Info<< " Lagrangian: " << lagrFileName << endl;
1246 
1247  wordList labelNames(sprayObjs.names(labelIOField::typeName));
1248  Info<< " labels :";
1249  print(Info, labelNames);
1250 
1251  wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
1252  Info<< " scalars :";
1253  print(Info, scalarNames);
1254 
1255  wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
1256  Info<< " vectors :";
1257  print(Info, vectorNames);
1258 
1259  // wordList sphereNames
1260  //(
1261  // sprayObjs.names
1262  // (
1263  // sphericalTensorIOField::typeName
1264  // )
1265  //);
1266  // Info<< " spherical tensors :";
1267  // print(Info, sphereNames);
1268  //
1269  // wordList symmNames
1270  //(
1271  // sprayObjs.names
1272  // (
1273  // symmTensorIOField::typeName
1274  // )
1275  //);
1276  // Info<< " symm tensors :";
1277  // print(Info, symmNames);
1278  //
1279  // wordList tensorNames
1280  //(
1281  // sprayObjs.names(tensorIOField::typeName)
1282  //);
1283  // Info<< " tensors :";
1284  // print(Info, tensorNames);
1285 
1286 
1287  // Load cloud positions
1288  passiveParticleCloud parcels(mesh, cloudDirs[cloudI]);
1289 
1290  // Get positions as pointField
1291  pointField positions(parcels.size());
1292  label n = 0;
1293  forAllConstIter(passiveParticleCloud, parcels, elmnt)
1294  {
1295  positions[n++] = elmnt().position();
1296  }
1297 
1298 
1299  string allVarNames = string("X Y Z");
1300  DynamicList<INTEGER4> allVarLocation;
1301  allVarLocation.append(ValueLocation_Nodal);
1302  allVarLocation.append(ValueLocation_Nodal);
1303  allVarLocation.append(ValueLocation_Nodal);
1304 
1305  tecplotWriter::getTecplotNames<label>
1306  (
1307  labelNames,
1308  ValueLocation_Nodal,
1309  allVarNames,
1310  allVarLocation
1311  );
1312 
1313  tecplotWriter::getTecplotNames<scalar>
1314  (
1315  scalarNames,
1316  ValueLocation_Nodal,
1317  allVarNames,
1318  allVarLocation
1319  );
1320 
1321  tecplotWriter::getTecplotNames<vector>
1322  (
1323  vectorNames,
1324  ValueLocation_Nodal,
1325  allVarNames,
1326  allVarLocation
1327  );
1328 
1329 
1330  tecplotWriter writer(runTime);
1331 
1332  writer.writeInit
1333  (
1334  runTime.caseName(),
1335  allVarNames,
1336  lagrFileName,
1337  DataFileType_Full
1338  );
1339 
1340  writer.writeOrderedZone
1341  (
1342  cloudDirs[cloudI],
1343  strandID++, // strandID,
1344  parcels.size(),
1345  allVarLocation
1346  );
1347 
1348  // Write coordinates
1349  writer.writeField(positions.component(0)());
1350  writer.writeField(positions.component(1)());
1351  writer.writeField(positions.component(2)());
1352 
1353  // labelFields
1354  forAll(labelNames, i)
1355  {
1356  IOField<label> fld
1357  (
1358  IOobject
1359  (
1360  labelNames[i],
1361  mesh.time().name(),
1362  lagrangian::cloud::prefix/cloudDirs[cloudI],
1363  mesh,
1364  IOobject::MUST_READ,
1365  IOobject::NO_WRITE,
1366  false
1367  )
1368  );
1369 
1370  scalarField sfld(fld.size());
1371  forAll(fld, j)
1372  {
1373  sfld[j] = scalar(fld[j]);
1374  }
1375  writer.writeField(sfld);
1376  }
1377  // scalarFields
1378  forAll(scalarNames, i)
1379  {
1380  IOField<scalar> fld
1381  (
1382  IOobject
1383  (
1384  scalarNames[i],
1385  mesh.time().name(),
1386  lagrangian::cloud::prefix/cloudDirs[cloudI],
1387  mesh,
1388  IOobject::MUST_READ,
1389  IOobject::NO_WRITE,
1390  false
1391  )
1392  );
1393  writer.writeField(fld);
1394  }
1395  // vectorFields
1396  forAll(vectorNames, i)
1397  {
1398  IOField<vector> fld
1399  (
1400  IOobject
1401  (
1402  vectorNames[i],
1403  mesh.time().name(),
1404  lagrangian::cloud::prefix/cloudDirs[cloudI],
1405  mesh,
1406  IOobject::MUST_READ,
1407  IOobject::NO_WRITE,
1408  false
1409  )
1410  );
1411  writer.writeField(fld);
1412  }
1413 
1414  writer.writeEnd();
1415  }
1416  }
1417  }
1418 
1419  Info<< "End\n" << endl;
1420 
1421  return 0;
1422 }
1423 
1424 
1425 // ************************************************************************* //
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:486
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:255
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:120
const word & name() const
Return const reference to name.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:443
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
int main(int argc, char *argv[])
Definition: financialFoam.C:44
volScalarField scalarField(fieldObject, mesh)
label patchi
static instantList timeDirs
Definition: globalFoam.H:44
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
List< word > wordList
A List of words.
Definition: fileName.H:54
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
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.
Definition: ReadFields.C:244
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
tmp< SurfaceField< Type > > linearInterpolate(const VolField< Type > &vf)
Definition: linear.H:108
List< instant > instantList
List of instants.
Definition: instantList.H:42
bool rmDir(const fileName &)
Remove a directory and its contents.
Definition: POSIX.C:1047
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
static const char nl
Definition: Ostream.H:267
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.
Definition: POSIX.C:662
objects
Foam::word regionName
Definition: setRegionName.H:1
Foam::argList args(argc, argv)
Operations on lists of strings.
mkDir(pdfPath)