FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

foamToVTK.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Application
00025     foamToVTK
00026 
00027 Description
00028     Legacy VTK file format writer.
00029 
00030     - handles volScalar, volVector, pointScalar, pointVector, surfaceScalar
00031       fields.
00032     - mesh topo changes.
00033     - both ascii and binary.
00034     - single time step writing.
00035     - write subset only.
00036     - automatic decomposition of cells; polygons on boundary undecomposed since
00037       handled by vtk.
00038 
00039 Usage
00040 
00041     - foamToVTK [OPTION]
00042 
00043 
00044     @param -ascii \n
00045     Write VTK data in ASCII format instead of binary.
00046 
00047     @param -mesh <name>\n
00048     Use a different mesh name (instead of -region)
00049 
00050     @param -fields <fields>\n
00051     Convert selected fields only. For example,
00052     @verbatim
00053          -fields "( p T U )"
00054     @endverbatim
00055     The quoting is required to avoid shell expansions and to pass the
00056     information as a single argument.
00057 
00058     @param -surfaceFields \n
00059     Write surfaceScalarFields (e.g., phi)
00060 
00061     @param -cellSet <name>\n
00062     @param -faceSet <name>\n
00063     @param -pointSet <name>\n.
00064     Restrict conversion to the cellSet, faceSet or pointSet.
00065 
00066     @param -nearCellValue \n
00067     Output cell value on patches instead of patch value itself.
00068 
00069     @param -noInternal \n
00070     Do not generate file for mesh, only for patches.
00071 
00072     @param -noPointValues \n
00073     No pointFields.
00074 
00075     @param -noFaceZones \n
00076     No faceZones.
00077 
00078     @param -noLinks \n
00079     (in parallel) do not link processor files to master.
00080 
00081     @param -allPatches \n
00082     Combine all patches into a single file.
00083 
00084     @param -excludePatches <patchNames>\n
00085     Specify patches to exclude. For example,
00086     @verbatim
00087          -excludePatches "( inlet_1 inlet_2 )"
00088     @endverbatim
00089     The quoting is required to avoid shell expansions and to pass the
00090     information as a single argument.
00091 
00092     @param -useTimeName \n
00093     use the time index in the VTK file name instead of the time index
00094 
00095     @param -case <dir> \n
00096     Case directory.
00097 
00098     @param -noZero \n
00099     Ignore time step 0.
00100 
00101     @param -constant \n
00102     Include the constant directory.
00103 
00104     @param -latestTime \n
00105     Only apply to the latest time step.
00106 
00107     @param -time <time>\n
00108     Apply only to specified time.
00109 
00110     @param -parallel \n
00111     Run in parallel.
00112 
00113     @param -help \n
00114     Display help message.
00115 
00116     @param -doc \n
00117     Display Doxygen API documentation page for this application.
00118 
00119     @param -srcDoc \n
00120     Display Doxygen source documentation page for this application.
00121 
00122 Note
00123     mesh subset is handled by vtkMesh. Slight inconsistency in
00124     interpolation: on the internal field it interpolates the whole volfield
00125     to the whole-mesh pointField and then selects only those values it
00126     needs for the subMesh (using the fvMeshSubset cellMap(), pointMap()
00127     functions). For the patches however it uses the
00128     fvMeshSubset.interpolate function to directly interpolate the
00129     whole-mesh values onto the subset patch.
00130 
00131 \*---------------------------------------------------------------------------*/
00132 
00133 #include <finiteVolume/fvCFD.H>
00134 #include <OpenFOAM/pointMesh.H>
00135 #include <finiteVolume/volPointInterpolation.H>
00136 #include <OpenFOAM/emptyPolyPatch.H>
00137 #include <OpenFOAM/labelIOField.H>
00138 #include <OpenFOAM/scalarIOField.H>
00139 #include <OpenFOAM/sphericalTensorIOField.H>
00140 #include <OpenFOAM/symmTensorIOField.H>
00141 #include <OpenFOAM/tensorIOField.H>
00142 #include <OpenFOAM/faceZoneMesh.H>
00143 #include <lagrangian/Cloud.H>
00144 #include <lagrangian/passiveParticle.H>
00145 
00146 #include "vtkMesh.H"
00147 #include "readFields.H"
00148 #include "writeFuns.H"
00149 
00150 #include "internalWriter.H"
00151 #include "patchWriter.H"
00152 #include "lagrangianWriter.H"
00153 
00154 #include "writeFaceSet.H"
00155 #include "writePointSet.H"
00156 #include "writePatchGeom.H"
00157 #include "writeSurfFields.H"
00158 
00159 
00160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00161 
00162 
00163 static const label VTK_TETRA      = 10;
00164 static const label VTK_PYRAMID    = 14;
00165 static const label VTK_WEDGE      = 13;
00166 static const label VTK_HEXAHEDRON = 12;
00167 
00168 
00169 template<class GeoField>
00170 void print(const char* msg, Ostream& os, const PtrList<GeoField>& flds)
00171 {
00172     if (flds.size())
00173     {
00174         os << msg;
00175         forAll(flds, i)
00176         {
00177             os<< ' ' << flds[i].name();
00178         }
00179         os << endl;
00180     }
00181 }
00182 
00183 
00184 void print(Ostream& os, const wordList& flds)
00185 {
00186     forAll(flds, i)
00187     {
00188         os<< ' ' << flds[i];
00189     }
00190     os << endl;
00191 }
00192 
00193 
00194 labelList getSelectedPatches
00195 (
00196     const polyBoundaryMesh& patches,
00197     const HashSet<word>& excludePatches
00198 )
00199 {
00200     DynamicList<label> patchIDs(patches.size());
00201 
00202     Info<< "Combining patches:" << endl;
00203 
00204     forAll(patches, patchI)
00205     {
00206         const polyPatch& pp = patches[patchI];
00207 
00208         if
00209         (
00210             isA<emptyPolyPatch>(pp)
00211             || (Pstream::parRun() && isA<processorPolyPatch>(pp))
00212         )
00213         {
00214             Info<< "    discarding empty/processor patch " << patchI
00215                 << " " << pp.name() << endl;
00216         }
00217         else if (!excludePatches.found(pp.name()))
00218         {
00219             patchIDs.append(patchI);
00220             Info<< "    patch " << patchI << " " << pp.name() << endl;
00221         }
00222     }
00223     return patchIDs.shrink();
00224 }
00225 
00226 
00227 
00228 
00229 // Main program:
00230 
00231 int main(int argc, char *argv[])
00232 {
00233     timeSelector::addOptions();
00234 
00235 #   include <OpenFOAM/addRegionOption.H>
00236 
00237     argList::validOptions.insert("fields", "fields");
00238     argList::validOptions.insert("cellSet", "cellSet name");
00239     argList::validOptions.insert("faceSet", "faceSet name");
00240     argList::validOptions.insert("pointSet", "pointSet name");
00241     argList::validOptions.insert("ascii","");
00242     argList::validOptions.insert("surfaceFields","");
00243     argList::validOptions.insert("nearCellValue","");
00244     argList::validOptions.insert("noInternal","");
00245     argList::validOptions.insert("noPointValues","");
00246     argList::validOptions.insert("allPatches","");
00247     argList::validOptions.insert("excludePatches","patches to exclude");
00248     argList::validOptions.insert("noFaceZones","");
00249     argList::validOptions.insert("noLinks","");
00250     argList::validOptions.insert("useTimeName","");
00251 
00252 #   include <OpenFOAM/setRootCase.H>
00253 #   include <OpenFOAM/createTime.H>
00254 
00255     bool doWriteInternal = !args.optionFound("noInternal");
00256     bool doFaceZones     = !args.optionFound("noFaceZones");
00257     bool doLinks         = !args.optionFound("noLinks");
00258     bool binary          = !args.optionFound("ascii");
00259     bool useTimeName     = args.optionFound("useTimeName");
00260 
00261     if (binary && (sizeof(floatScalar) != 4 || sizeof(label) != 4))
00262     {
00263         FatalErrorIn(args.executable())
00264             << "floatScalar and/or label are not 4 bytes in size" << nl
00265             << "Hence cannot use binary VTK format. Please use -ascii"
00266             << exit(FatalError);
00267     }
00268 
00269     bool nearCellValue = args.optionFound("nearCellValue");
00270 
00271     if (nearCellValue)
00272     {
00273         WarningIn(args.executable())
00274             << "Using neighbouring cell value instead of patch value"
00275             << nl << endl;
00276     }
00277 
00278     bool noPointValues = args.optionFound("noPointValues");
00279 
00280     if (noPointValues)
00281     {
00282         WarningIn(args.executable())
00283             << "Outputting cell values only" << nl << endl;
00284     }
00285 
00286     bool allPatches = args.optionFound("allPatches");
00287 
00288     HashSet<word> excludePatches;
00289     if (args.optionFound("excludePatches"))
00290     {
00291         args.optionLookup("excludePatches")() >> excludePatches;
00292 
00293         Info<< "Not including patches " << excludePatches << nl << endl;
00294     }
00295 
00296     word cellSetName;
00297     string vtkName;
00298 
00299     if (args.optionFound("cellSet"))
00300     {
00301         cellSetName = args.option("cellSet");
00302         vtkName = cellSetName;
00303     }
00304     else if (Pstream::parRun())
00305     {
00306         // Strip off leading casename, leaving just processor_DDD ending.
00307         vtkName = runTime.caseName();
00308 
00309         string::size_type i = vtkName.rfind("processor");
00310 
00311         if (i != string::npos)
00312         {
00313             vtkName = vtkName.substr(i);
00314         }
00315     }
00316     else
00317     {
00318         vtkName = runTime.caseName();
00319     }
00320 
00321 
00322     instantList timeDirs = timeSelector::select0(runTime, args);
00323 
00324 #   include <OpenFOAM/createNamedMesh.H>
00325 
00326     // VTK/ directory in the case
00327     fileName fvPath(runTime.path()/"VTK");
00328     // Directory of mesh (region0 gets filtered out)
00329     fileName regionPrefix = "";
00330 
00331     if (regionName != polyMesh::defaultRegion)
00332     {
00333         fvPath = fvPath/regionName;
00334         regionPrefix = regionName;
00335     }
00336 
00337     if (isDir(fvPath))
00338     {
00339         if
00340         (
00341             args.optionFound("time")
00342          || args.optionFound("latestTime")
00343          || cellSetName.size()
00344          || regionName != polyMesh::defaultRegion
00345         )
00346         {
00347             Info<< "Keeping old VTK files in " << fvPath << nl << endl;
00348         }
00349         else
00350         {
00351             Info<< "Deleting old VTK files in " << fvPath << nl << endl;
00352 
00353             rmDir(fvPath);
00354         }
00355     }
00356 
00357     mkDir(fvPath);
00358 
00359 
00360     // mesh wrapper; does subsetting and decomposition
00361     vtkMesh vMesh(mesh, cellSetName);
00362 
00363 
00364     // Scan for all possible lagrangian clouds
00365     HashSet<fileName> allCloudDirs;
00366     forAll(timeDirs, timeI)
00367     {
00368         runTime.setTime(timeDirs[timeI], timeI);
00369         fileNameList cloudDirs
00370         (
00371             readDir
00372             (
00373                 runTime.timePath()/regionPrefix/cloud::prefix,
00374                 fileName::DIRECTORY
00375             )
00376         );
00377         forAll(cloudDirs, i)
00378         {
00379             IOobjectList sprayObjs
00380             (
00381                 mesh,
00382                 runTime.timeName(),
00383                 cloud::prefix/cloudDirs[i]
00384             );
00385 
00386             IOobject* positionsPtr = sprayObjs.lookup("positions");
00387 
00388             if (positionsPtr)
00389             {
00390                 if (allCloudDirs.insert(cloudDirs[i]))
00391                 {
00392                     Info<< "At time: " << runTime.timeName()
00393                         << " detected cloud directory : " << cloudDirs[i]
00394                         << endl;
00395                 }
00396             }
00397         }
00398     }
00399 
00400 
00401     forAll(timeDirs, timeI)
00402     {
00403         runTime.setTime(timeDirs[timeI], timeI);
00404 
00405         Info<< "Time: " << runTime.timeName() << endl;
00406 
00407         word timeDesc = "";
00408         if (useTimeName)
00409         {
00410             timeDesc = runTime.timeName();
00411         }
00412         else
00413         {
00414             timeDesc = name(runTime.timeIndex());
00415         }
00416 
00417         // Check for new polyMesh/ and update mesh, fvMeshSubset and cell
00418         // decomposition.
00419         polyMesh::readUpdateState meshState = vMesh.readUpdate();
00420 
00421         const fvMesh& mesh = vMesh.mesh();
00422 
00423         if
00424         (
00425             meshState == polyMesh::TOPO_CHANGE
00426          || meshState == polyMesh::TOPO_PATCH_CHANGE
00427         )
00428         {
00429             Info<< "    Read new mesh" << nl << endl;
00430         }
00431 
00432 
00433         // If faceSet: write faceSet only (as polydata)
00434         if (args.optionFound("faceSet"))
00435         {
00436             // Load the faceSet
00437             faceSet set(mesh, args.option("faceSet"));
00438 
00439             // Filename as if patch with same name.
00440             mkDir(fvPath/set.name());
00441 
00442             fileName patchFileName
00443             (
00444                 fvPath/set.name()/set.name()
00445               + "_"
00446               + timeDesc
00447               + ".vtk"
00448             );
00449 
00450             Info<< "    FaceSet   : " << patchFileName << endl;
00451 
00452             writeFaceSet(binary, vMesh, set, patchFileName);
00453 
00454             continue;
00455         }
00456         // If pointSet: write pointSet only (as polydata)
00457         if (args.optionFound("pointSet"))
00458         {
00459             // Load the pointSet
00460             pointSet set(mesh, args.option("pointSet"));
00461 
00462             // Filename as if patch with same name.
00463             mkDir(fvPath/set.name());
00464 
00465             fileName patchFileName
00466             (
00467                 fvPath/set.name()/set.name()
00468               + "_"
00469               + timeDesc
00470               + ".vtk"
00471             );
00472 
00473             Info<< "    pointSet   : " << patchFileName << endl;
00474 
00475             writePointSet(binary, vMesh, set, patchFileName);
00476 
00477             continue;
00478         }
00479 
00480 
00481         // Search for list of objects for this time
00482         IOobjectList objects(mesh, runTime.timeName());
00483 
00484         HashSet<word> selectedFields;
00485         if (args.optionFound("fields"))
00486         {
00487             args.optionLookup("fields")() >> selectedFields;
00488         }
00489 
00490         // Construct the vol fields (on the original mesh if subsetted)
00491 
00492         PtrList<volScalarField> vsf;
00493         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vsf);
00494         print("    volScalarFields            :", Info, vsf);
00495 
00496         PtrList<volVectorField> vvf;
00497         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vvf);
00498         print("    volVectorFields            :", Info, vvf);
00499 
00500         PtrList<volSphericalTensorField> vSpheretf;
00501         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSpheretf);
00502         print("    volSphericalTensorFields   :", Info, vSpheretf);
00503 
00504         PtrList<volSymmTensorField> vSymmtf;
00505         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vSymmtf);
00506         print("    volSymmTensorFields        :", Info, vSymmtf);
00507 
00508         PtrList<volTensorField> vtf;
00509         readFields(vMesh, vMesh.baseMesh(), objects, selectedFields, vtf);
00510         print("    volTensorFields            :", Info, vtf);
00511 
00512         label nVolFields =
00513                 vsf.size()
00514               + vvf.size()
00515               + vSpheretf.size()
00516               + vSymmtf.size()
00517               + vtf.size();
00518 
00519 
00520         // Construct pointMesh only if nessecary since constructs edge
00521         // addressing (expensive on polyhedral meshes)
00522         if (noPointValues)
00523         {
00524             Info<< "    pointScalarFields : switched off"
00525                 << " (\"-noPointValues\" option)\n";
00526             Info<< "    pointVectorFields : switched off"
00527                 << " (\"-noPointValues\" option)\n";
00528         }
00529 
00530         PtrList<pointScalarField> psf;
00531         PtrList<pointVectorField> pvf;
00532         PtrList<pointSphericalTensorField> pSpheretf;
00533         PtrList<pointSymmTensorField> pSymmtf;
00534         PtrList<pointTensorField> ptf;
00535 
00536         if (!noPointValues)
00537         {
00538             readFields
00539             (
00540                 vMesh,
00541                 pointMesh::New(vMesh.baseMesh()),
00542                 objects,
00543                 selectedFields,
00544                 psf
00545             );
00546             print("    pointScalarFields          :", Info, psf);
00547 
00548             readFields
00549             (
00550                 vMesh,
00551                 pointMesh::New(vMesh.baseMesh()),
00552                 objects,
00553                 selectedFields,
00554                 pvf
00555             );
00556             print("    pointVectorFields          :", Info, pvf);
00557 
00558             readFields
00559             (
00560                 vMesh,
00561                 pointMesh::New(vMesh.baseMesh()),
00562                 objects,
00563                 selectedFields,
00564                 pSpheretf
00565             );
00566             print("    pointSphericalTensorFields :", Info, pSpheretf);
00567 
00568             readFields
00569             (
00570                 vMesh,
00571                 pointMesh::New(vMesh.baseMesh()),
00572                 objects,
00573                 selectedFields,
00574                 pSymmtf
00575             );
00576             print("    pointSymmTensorFields      :", Info, pSymmtf);
00577 
00578             readFields
00579             (
00580                 vMesh,
00581                 pointMesh::New(vMesh.baseMesh()),
00582                 objects,
00583                 selectedFields,
00584                 ptf
00585             );
00586             print("    pointTensorFields          :", Info, ptf);
00587         }
00588         Info<< endl;
00589 
00590         label nPointFields =
00591             psf.size()
00592           + pvf.size()
00593           + pSpheretf.size()
00594           + pSymmtf.size()
00595           + ptf.size();
00596 
00597         if (doWriteInternal)
00598         {
00599             //
00600             // Create file and write header
00601             //
00602             fileName vtkFileName
00603             (
00604                 fvPath/vtkName
00605               + "_"
00606               + timeDesc
00607               + ".vtk"
00608             );
00609 
00610             Info<< "    Internal  : " << vtkFileName << endl;
00611 
00612             // Write mesh
00613             internalWriter writer(vMesh, binary, vtkFileName);
00614 
00615             // VolFields + cellID
00616             writeFuns::writeCellDataHeader
00617             (
00618                 writer.os(),
00619                 vMesh.nFieldCells(),
00620                 1+nVolFields
00621             );
00622 
00623             // Write cellID field
00624             writer.writeCellIDs();
00625 
00626             // Write volFields
00627             writer.write(vsf);
00628             writer.write(vvf);
00629             writer.write(vSpheretf);
00630             writer.write(vSymmtf);
00631             writer.write(vtf);
00632 
00633             if (!noPointValues)
00634             {
00635                 writeFuns::writePointDataHeader
00636                 (
00637                     writer.os(),
00638                     vMesh.nFieldPoints(),
00639                     nVolFields+nPointFields
00640                 );
00641 
00642                 // pointFields
00643                 writer.write(psf);
00644                 writer.write(pvf);
00645                 writer.write(pSpheretf);
00646                 writer.write(pSymmtf);
00647                 writer.write(ptf);
00648 
00649                 // Interpolated volFields
00650                 volPointInterpolation pInterp(mesh);
00651                 writer.write(pInterp, vsf);
00652                 writer.write(pInterp, vvf);
00653                 writer.write(pInterp, vSpheretf);
00654                 writer.write(pInterp, vSymmtf);
00655                 writer.write(pInterp, vtf);
00656             }
00657         }
00658 
00659         //---------------------------------------------------------------------
00660         //
00661         // Write surface fields
00662         //
00663         //---------------------------------------------------------------------
00664 
00665         if (args.optionFound("surfaceFields"))
00666         {
00667             PtrList<surfaceScalarField> ssf;
00668             readFields
00669             (
00670                 vMesh,
00671                 vMesh.baseMesh(),
00672                 objects,
00673                 selectedFields,
00674                 ssf
00675             );
00676             print("    surfScalarFields  :", Info, ssf);
00677 
00678             PtrList<surfaceVectorField> svf;
00679             readFields
00680             (
00681                 vMesh,
00682                 vMesh.baseMesh(),
00683                 objects,
00684                 selectedFields,
00685                 svf
00686             );
00687             print("    surfVectorFields  :", Info, svf);
00688 
00689             if (ssf.size() + svf.size() > 0)
00690             {
00691                 // Rework the scalar fields into vectorfields.
00692                 label sz = svf.size();
00693 
00694                 svf.setSize(sz+ssf.size());
00695 
00696                 surfaceVectorField n(mesh.Sf()/mesh.magSf());
00697 
00698                 forAll(ssf, i)
00699                 {
00700                     svf.set(sz+i, ssf[i]*n);
00701                     svf[sz+i].rename(ssf[i].name());
00702                 }
00703                 ssf.clear();
00704 
00705                 mkDir(fvPath / "surfaceFields");
00706 
00707                 fileName surfFileName
00708                 (
00709                     fvPath
00710                    /"surfaceFields"
00711                    /"surfaceFields"
00712                    + "_"
00713                    + timeDesc
00714                    + ".vtk"
00715                 );
00716 
00717                 writeSurfFields
00718                 (
00719                     binary,
00720                     vMesh,
00721                     surfFileName,
00722                     svf
00723                 );
00724             }
00725         }
00726 
00727 
00728         //---------------------------------------------------------------------
00729         //
00730         // Write patches (POLYDATA file, one for each patch)
00731         //
00732         //---------------------------------------------------------------------
00733 
00734         const polyBoundaryMesh& patches = mesh.boundaryMesh();
00735 
00736         if (allPatches)
00737         {
00738             mkDir(fvPath/"allPatches");
00739 
00740             fileName patchFileName;
00741 
00742             if (vMesh.useSubMesh())
00743             {
00744                 patchFileName =
00745                     fvPath/"allPatches"/cellSetName
00746                   + "_"
00747                   + timeDesc
00748                   + ".vtk";
00749             }
00750             else
00751             {
00752                 patchFileName =
00753                     fvPath/"allPatches"/"allPatches"
00754                   + "_"
00755                   + timeDesc
00756                   + ".vtk";
00757             }
00758 
00759             Info<< "    Combined patches     : " << patchFileName << endl;
00760 
00761             patchWriter writer
00762             (
00763                 vMesh,
00764                 binary,
00765                 nearCellValue,
00766                 patchFileName,
00767                 getSelectedPatches(patches, excludePatches)
00768             );
00769 
00770             // VolFields + patchID
00771             writeFuns::writeCellDataHeader
00772             (
00773                 writer.os(),
00774                 writer.nFaces(),
00775                 1+nVolFields
00776             );
00777 
00778             // Write patchID field
00779             writer.writePatchIDs();
00780 
00781             // Write volFields
00782             writer.write(vsf);
00783             writer.write(vvf);
00784             writer.write(vSpheretf);
00785             writer.write(vSymmtf);
00786             writer.write(vtf);
00787 
00788             if (!noPointValues)
00789             {
00790                 writeFuns::writePointDataHeader
00791                 (
00792                     writer.os(),
00793                     writer.nPoints(),
00794                     nPointFields
00795                 );
00796 
00797                 // Write pointFields
00798                 writer.write(psf);
00799                 writer.write(pvf);
00800                 writer.write(pSpheretf);
00801                 writer.write(pSymmtf);
00802                 writer.write(ptf);
00803 
00804                 // no interpolated volFields since I cannot be bothered to
00805                 // create the patchInterpolation for all subpatches.
00806             }
00807         }
00808         else
00809         {
00810             forAll(patches, patchI)
00811             {
00812                 const polyPatch& pp = patches[patchI];
00813 
00814                 if (!excludePatches.found(pp.name()))
00815                 {
00816                     mkDir(fvPath/pp.name());
00817 
00818                     fileName patchFileName;
00819 
00820                     if (vMesh.useSubMesh())
00821                     {
00822                         patchFileName =
00823                             fvPath/pp.name()/cellSetName
00824                           + "_"
00825                           + timeDesc
00826                           + ".vtk";
00827                     }
00828                     else
00829                     {
00830                         patchFileName =
00831                             fvPath/pp.name()/pp.name()
00832                           + "_"
00833                           + timeDesc
00834                           + ".vtk";
00835                     }
00836 
00837                     Info<< "    Patch     : " << patchFileName << endl;
00838 
00839                     patchWriter writer
00840                     (
00841                         vMesh,
00842                         binary,
00843                         nearCellValue,
00844                         patchFileName,
00845                         labelList(1, patchI)
00846                     );
00847 
00848                     if (!isA<emptyPolyPatch>(pp))
00849                     {
00850                         // VolFields + patchID
00851                         writeFuns::writeCellDataHeader
00852                         (
00853                             writer.os(),
00854                             writer.nFaces(),
00855                             1+nVolFields
00856                         );
00857 
00858                         // Write patchID field
00859                         writer.writePatchIDs();
00860 
00861                         // Write volFields
00862                         writer.write(vsf);
00863                         writer.write(vvf);
00864                         writer.write(vSpheretf);
00865                         writer.write(vSymmtf);
00866                         writer.write(vtf);
00867 
00868                         if (!noPointValues)
00869                         {
00870                             writeFuns::writePointDataHeader
00871                             (
00872                                 writer.os(),
00873                                 writer.nPoints(),
00874                                 nVolFields
00875                               + nPointFields
00876                             );
00877 
00878                             // Write pointFields
00879                             writer.write(psf);
00880                             writer.write(pvf);
00881                             writer.write(pSpheretf);
00882                             writer.write(pSymmtf);
00883                             writer.write(ptf);
00884 
00885                             PrimitivePatchInterpolation<primitivePatch> pInter
00886                             (
00887                                 pp
00888                             );
00889 
00890                             // Write interpolated volFields
00891                             writer.write(pInter, vsf);
00892                             writer.write(pInter, vvf);
00893                             writer.write(pInter, vSpheretf);
00894                             writer.write(pInter, vSymmtf);
00895                             writer.write(pInter, vtf);
00896                         }
00897                     }
00898                 }
00899             }
00900         }
00901 
00902         //---------------------------------------------------------------------
00903         //
00904         // Write faceZones (POLYDATA file, one for each zone)
00905         //
00906         //---------------------------------------------------------------------
00907 
00908         if (doFaceZones)
00909         {
00910             const faceZoneMesh& zones = mesh.faceZones();
00911 
00912             forAll(zones, zoneI)
00913             {
00914                 const faceZone& pp = zones[zoneI];
00915 
00916                 mkDir(fvPath/pp.name());
00917 
00918                 fileName patchFileName;
00919 
00920                 if (vMesh.useSubMesh())
00921                 {
00922                     patchFileName =
00923                         fvPath/pp.name()/cellSetName
00924                       + "_"
00925                       + timeDesc
00926                       + ".vtk";
00927                 }
00928                 else
00929                 {
00930                     patchFileName =
00931                         fvPath/pp.name()/pp.name()
00932                       + "_"
00933                       + timeDesc
00934                       + ".vtk";
00935                 }
00936 
00937                 Info<< "    FaceZone  : " << patchFileName << endl;
00938 
00939                 std::ofstream str(patchFileName.c_str());
00940 
00941                 writeFuns::writeHeader(str, binary, pp.name());
00942                 str << "DATASET POLYDATA" << std::endl;
00943 
00944                 writePatchGeom
00945                 (
00946                     binary,
00947                     pp().localFaces(),
00948                     pp().localPoints(),
00949                     str
00950                 );
00951             }
00952         }
00953 
00954 
00955 
00956         //---------------------------------------------------------------------
00957         //
00958         // Write lagrangian data
00959         //
00960         //---------------------------------------------------------------------
00961 
00962         forAllConstIter(HashSet<fileName>, allCloudDirs, iter)
00963         {
00964             const fileName& cloudName = iter.key();
00965 
00966             // Always create the cloud directory.
00967             mkDir(fvPath/cloud::prefix/cloudName);
00968 
00969             fileName lagrFileName
00970             (
00971                 fvPath/cloud::prefix/cloudName/cloudName
00972               + "_" + timeDesc + ".vtk"
00973             );
00974 
00975             Info<< "    Lagrangian: " << lagrFileName << endl;
00976 
00977 
00978             IOobjectList sprayObjs
00979             (
00980                 mesh,
00981                 runTime.timeName(),
00982                 cloud::prefix/cloudName
00983             );
00984 
00985             IOobject* positionsPtr = sprayObjs.lookup("positions");
00986 
00987             if (positionsPtr)
00988             {
00989                 wordList labelNames(sprayObjs.names(labelIOField::typeName));
00990                 Info<< "        labels            :";
00991                 print(Info, labelNames);
00992 
00993                 wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
00994                 Info<< "        scalars           :";
00995                 print(Info, scalarNames);
00996 
00997                 wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
00998                 Info<< "        vectors           :";
00999                 print(Info, vectorNames);
01000 
01001                 wordList sphereNames
01002                 (
01003                     sprayObjs.names
01004                     (
01005                         sphericalTensorIOField::typeName
01006                     )
01007                 );
01008                 Info<< "        spherical tensors :";
01009                 print(Info, sphereNames);
01010 
01011                 wordList symmNames
01012                 (
01013                     sprayObjs.names
01014                     (
01015                         symmTensorIOField::typeName
01016                     )
01017                 );
01018                 Info<< "        symm tensors      :";
01019                 print(Info, symmNames);
01020 
01021                 wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
01022                 Info<< "        tensors           :";
01023                 print(Info, tensorNames);
01024 
01025                 lagrangianWriter writer
01026                 (
01027                     vMesh,
01028                     binary,
01029                     lagrFileName,
01030                     cloudName,
01031                     false
01032                 );
01033 
01034                 // Write number of fields
01035                 writer.writeParcelHeader
01036                 (
01037                     labelNames.size()
01038                   + scalarNames.size()
01039                   + vectorNames.size()
01040                   + sphereNames.size()
01041                   + symmNames.size()
01042                   + tensorNames.size()
01043                 );
01044 
01045                 // Fields
01046                 writer.writeIOField<label>(labelNames);
01047                 writer.writeIOField<scalar>(scalarNames);
01048                 writer.writeIOField<vector>(vectorNames);
01049                 writer.writeIOField<sphericalTensor>(sphereNames);
01050                 writer.writeIOField<symmTensor>(symmNames);
01051                 writer.writeIOField<tensor>(tensorNames);
01052             }
01053             else
01054             {
01055                 lagrangianWriter writer
01056                 (
01057                     vMesh,
01058                     binary,
01059                     lagrFileName,
01060                     cloudName,
01061                     true
01062                 );
01063 
01064                 // Write number of fields
01065                 writer.writeParcelHeader(0);
01066             }
01067         }
01068     }
01069 
01070 
01071     //---------------------------------------------------------------------
01072     //
01073     // Link parallel outputs back to undecomposed case for ease of loading
01074     //
01075     //---------------------------------------------------------------------
01076 
01077     if (Pstream::parRun() && doLinks)
01078     {
01079         mkDir(runTime.path()/".."/"VTK");
01080         chDir(runTime.path()/".."/"VTK");
01081 
01082         Info<< "Linking all processor files to " << runTime.path()/".."/"VTK"
01083             << endl;
01084 
01085         // Get list of vtk files
01086         fileName procVTK
01087         (
01088             fileName("..")
01089            /"processor" + name(Pstream::myProcNo())
01090            /"VTK"
01091         );
01092 
01093         fileNameList dirs(readDir(procVTK, fileName::DIRECTORY));
01094         label sz = dirs.size();
01095         dirs.setSize(sz+1);
01096         dirs[sz] = ".";
01097 
01098         forAll(dirs, i)
01099         {
01100             fileNameList subFiles(readDir(procVTK/dirs[i], fileName::FILE));
01101 
01102             forAll(subFiles, j)
01103             {
01104                 fileName procFile(procVTK/dirs[i]/subFiles[j]);
01105 
01106                 if (exists(procFile))
01107                 {
01108                     string cmd
01109                     (
01110                         "ln -s "
01111                       + procFile
01112                       + " "
01113                       + "processor"
01114                       + name(Pstream::myProcNo())
01115                       + "_"
01116                       + procFile.name()
01117                     );
01118                     if (system(cmd.c_str()) == -1)
01119                     {
01120                         WarningIn(args.executable())
01121                             << "Could not execute command " << cmd << endl;
01122                     }
01123                 }
01124             }
01125         }
01126     }
01127 
01128     Info<< "End\n" << endl;
01129 
01130     return 0;
01131 }
01132 
01133 
01134 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines