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