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