2010-05-12 13:27:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
2013-12-11 16:09:41 +00:00
|
|
|
\\ / F ield | foam-extend: Open Source CFD
|
2018-05-29 07:35:20 +00:00
|
|
|
\\ / O peration | Version: 4.1
|
2015-05-17 13:32:07 +00:00
|
|
|
\\ / A nd | Web: http://www.foam-extend.org
|
|
|
|
\\/ M anipulation | For copyright notice see file Copyright
|
2010-05-12 13:27:55 +00:00
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
2013-12-11 16:09:41 +00:00
|
|
|
This file is part of foam-extend.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is free software: you can redistribute it and/or modify it
|
2010-05-12 13:27:55 +00:00
|
|
|
under the terms of the GNU General Public License as published by the
|
2013-12-11 16:09:41 +00:00
|
|
|
Free Software Foundation, either version 3 of the License, or (at your
|
2010-05-12 13:27:55 +00:00
|
|
|
option) any later version.
|
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is distributed in the hope that it will be useful, but
|
|
|
|
WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
|
|
|
|
General Public License for more details.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
2013-12-11 16:09:41 +00:00
|
|
|
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
Application
|
|
|
|
foamToVTK
|
|
|
|
|
|
|
|
Description
|
|
|
|
Legacy VTK file format writer.
|
|
|
|
|
|
|
|
- handles volScalar, volVector, pointScalar, pointVector, surfaceScalar
|
|
|
|
fields.
|
|
|
|
- mesh topo changes.
|
|
|
|
- both ascii and binary.
|
|
|
|
- single time step writing.
|
|
|
|
- write subset only.
|
|
|
|
- automatic decomposition of cells; polygons on boundary undecomposed since
|
|
|
|
handled by vtk.
|
|
|
|
|
|
|
|
Usage
|
|
|
|
|
|
|
|
- foamToVTK [OPTION]
|
|
|
|
|
|
|
|
|
|
|
|
@param -ascii \n
|
|
|
|
Write VTK data in ASCII format instead of binary.
|
|
|
|
|
|
|
|
@param -mesh \<name\>\n
|
|
|
|
Use a different mesh name (instead of -region)
|
|
|
|
|
|
|
|
@param -fields \<fields\>\n
|
|
|
|
Convert selected fields only. For example,
|
|
|
|
@verbatim
|
|
|
|
-fields "( p T U )"
|
|
|
|
@endverbatim
|
|
|
|
The quoting is required to avoid shell expansions and to pass the
|
|
|
|
information as a single argument.
|
|
|
|
|
|
|
|
@param -surfaceFields \n
|
|
|
|
Write surfaceScalarFields (e.g., phi)
|
|
|
|
|
|
|
|
@param -cellSet \<name\>\n
|
|
|
|
@param -faceSet \<name\>\n
|
|
|
|
@param -pointSet \<name\>\n
|
|
|
|
Restrict conversion to the cellSet, faceSet or pointSet.
|
|
|
|
|
|
|
|
@param -nearCellValue \n
|
|
|
|
Output cell value on patches instead of patch value itself
|
|
|
|
|
|
|
|
@param -noInternal \n
|
|
|
|
Do not generate file for mesh, only for patches
|
|
|
|
|
|
|
|
@param -noPointValues \n
|
|
|
|
No pointFields
|
|
|
|
|
|
|
|
@param -noFaceZones \n
|
|
|
|
No faceZones
|
|
|
|
|
|
|
|
@param -noLinks \n
|
|
|
|
(in parallel) do not link processor files to master
|
|
|
|
|
|
|
|
@param -allPatches \n
|
|
|
|
Combine all patches into a single file
|
|
|
|
|
2011-02-25 15:25:43 +00:00
|
|
|
@param -allWallPatches \n
|
|
|
|
Combine all wall patches into a single file
|
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
@param -excludePatches \<patchNames\>\n
|
|
|
|
Specify patches to exclude. For example,
|
|
|
|
@verbatim
|
|
|
|
-excludePatches "( inlet_1 inlet_2 )"
|
|
|
|
@endverbatim
|
|
|
|
The quoting is required to avoid shell expansions and to pass the
|
|
|
|
information as a single argument.
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
@param -useTimeName \n
|
|
|
|
use the time index in the VTK file name instead of the time index
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
Note
|
|
|
|
mesh subset is handled by vtkMesh. Slight inconsistency in
|
|
|
|
interpolation: on the internal field it interpolates the whole volfield
|
|
|
|
to the whole-mesh pointField and then selects only those values it
|
|
|
|
needs for the subMesh (using the fvMeshSubset cellMap(), pointMap()
|
|
|
|
functions). For the patches however it uses the
|
|
|
|
fvMeshSubset.interpolate function to directly interpolate the
|
|
|
|
whole-mesh values onto the subset patch.
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
|
|
|
#include "pointMesh.H"
|
|
|
|
#include "volPointInterpolation.H"
|
|
|
|
#include "emptyPolyPatch.H"
|
|
|
|
#include "labelIOField.H"
|
|
|
|
#include "scalarIOField.H"
|
|
|
|
#include "sphericalTensorIOField.H"
|
|
|
|
#include "symmTensorIOField.H"
|
|
|
|
#include "tensorIOField.H"
|
|
|
|
#include "faceZoneMesh.H"
|
2015-06-29 17:34:30 +00:00
|
|
|
#include "CloudTemplate.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
#include "passiveParticle.H"
|
|
|
|
|
|
|
|
#include "faCFD.H"
|
|
|
|
|
|
|
|
#include "vtkMesh.H"
|
|
|
|
#include "readFields.H"
|
|
|
|
#include "writeFuns.H"
|
|
|
|
|
|
|
|
#include "internalWriter.H"
|
|
|
|
#include "patchWriter.H"
|
|
|
|
#include "faMeshWriter.H"
|
|
|
|
#include "lagrangianWriter.H"
|
|
|
|
|
|
|
|
#include "writeFaceSet.H"
|
|
|
|
#include "writePointSet.H"
|
|
|
|
#include "writePatchGeom.H"
|
|
|
|
#include "writeSurfFields.H"
|
|
|
|
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
|
2011-02-25 15:25:43 +00:00
|
|
|
static const label VTK_TETRA = 10;
|
|
|
|
static const label VTK_PYRAMID = 14;
|
|
|
|
static const label VTK_WEDGE = 13;
|
2010-05-12 13:27:55 +00:00
|
|
|
static const label VTK_HEXAHEDRON = 12;
|
|
|
|
|
|
|
|
|
|
|
|
template<class GeoField>
|
|
|
|
void print(const char* msg, Ostream& os, const PtrList<GeoField>& flds)
|
|
|
|
{
|
|
|
|
if (flds.size())
|
|
|
|
{
|
|
|
|
os << msg;
|
|
|
|
forAll(flds, i)
|
|
|
|
{
|
|
|
|
os<< ' ' << flds[i].name();
|
|
|
|
}
|
|
|
|
os << endl;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
void print(Ostream& os, const wordList& flds)
|
|
|
|
{
|
|
|
|
forAll(flds, i)
|
|
|
|
{
|
|
|
|
os<< ' ' << flds[i];
|
|
|
|
}
|
|
|
|
os << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
labelList getSelectedPatches
|
|
|
|
(
|
|
|
|
const polyBoundaryMesh& patches,
|
|
|
|
const HashSet<word>& excludePatches
|
|
|
|
)
|
|
|
|
{
|
|
|
|
DynamicList<label> patchIDs(patches.size());
|
|
|
|
|
|
|
|
Info<< "Combining patches:" << endl;
|
|
|
|
|
|
|
|
forAll(patches, patchI)
|
|
|
|
{
|
|
|
|
const polyPatch& pp = patches[patchI];
|
|
|
|
|
|
|
|
if
|
|
|
|
(
|
2010-08-26 14:22:03 +00:00
|
|
|
isA<emptyPolyPatch>(pp)
|
|
|
|
|| (Pstream::parRun() && isA<processorPolyPatch>(pp))
|
2010-05-12 13:27:55 +00:00
|
|
|
)
|
|
|
|
{
|
|
|
|
Info<< " discarding empty/processor patch " << patchI
|
|
|
|
<< " " << pp.name() << endl;
|
|
|
|
}
|
|
|
|
else if (!excludePatches.found(pp.name()))
|
|
|
|
{
|
|
|
|
patchIDs.append(patchI);
|
|
|
|
Info<< " patch " << patchI << " " << pp.name() << endl;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
return patchIDs.shrink();
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
// Main program:
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
timeSelector::addOptions();
|
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
# include "addRegionOption.H"
|
2010-08-26 14:22:03 +00:00
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
argList::validOptions.insert("fields", "fields");
|
|
|
|
argList::validOptions.insert("cellSet", "cellSet name");
|
|
|
|
argList::validOptions.insert("faceSet", "faceSet name");
|
|
|
|
argList::validOptions.insert("pointSet", "pointSet name");
|
|
|
|
argList::validOptions.insert("ascii","");
|
|
|
|
argList::validOptions.insert("surfaceFields","");
|
|
|
|
argList::validOptions.insert("nearCellValue","");
|
|
|
|
argList::validOptions.insert("noInternal","");
|
|
|
|
argList::validOptions.insert("noPointValues","");
|
|
|
|
argList::validOptions.insert("allPatches","");
|
2011-02-25 15:25:43 +00:00
|
|
|
argList::validOptions.insert("allWallPatches","");
|
2010-05-12 13:27:55 +00:00
|
|
|
argList::validOptions.insert("excludePatches","patches to exclude");
|
|
|
|
argList::validOptions.insert("noFaceZones","");
|
|
|
|
argList::validOptions.insert("faMesh","");
|
|
|
|
argList::validOptions.insert("noLinks","");
|
2010-08-26 14:22:03 +00:00
|
|
|
argList::validOptions.insert("useTimeName","");
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
# include "setRootCase.H"
|
|
|
|
# include "createTime.H"
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
bool doWriteInternal = !args.optionFound("noInternal");
|
|
|
|
bool doFaceZones = !args.optionFound("noFaceZones");
|
|
|
|
bool doLinks = !args.optionFound("noLinks");
|
|
|
|
bool binary = !args.optionFound("ascii");
|
|
|
|
bool useTimeName = args.optionFound("useTimeName");
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
if (binary && (sizeof(floatScalar) != 4 || sizeof(label) != 4))
|
|
|
|
{
|
|
|
|
FatalErrorIn(args.executable())
|
|
|
|
<< "floatScalar and/or label are not 4 bytes in size" << nl
|
|
|
|
<< "Hence cannot use binary VTK format. Please use -ascii"
|
|
|
|
<< exit(FatalError);
|
|
|
|
}
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
bool nearCellValue = args.optionFound("nearCellValue");
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
if (nearCellValue)
|
|
|
|
{
|
2011-02-25 15:25:43 +00:00
|
|
|
Info<< "Using neighbouring cell value instead of patch value"
|
2010-05-12 13:27:55 +00:00
|
|
|
<< nl << endl;
|
|
|
|
}
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
bool noPointValues = args.optionFound("noPointValues");
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
if (noPointValues)
|
|
|
|
{
|
2011-02-25 15:25:43 +00:00
|
|
|
Info<< "Outputting cell values only" << nl << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
bool allPatches = args.optionFound("allPatches");
|
2011-02-25 15:25:43 +00:00
|
|
|
bool allWallPatches = args.optionFound("allWallPatches");
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
HashSet<word> excludePatches;
|
2010-08-26 14:22:03 +00:00
|
|
|
if (args.optionFound("excludePatches"))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
args.optionLookup("excludePatches")() >> excludePatches;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
Info<< "Not including patches " << excludePatches << nl << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
word cellSetName;
|
|
|
|
string vtkName;
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
if (args.optionFound("cellSet"))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
cellSetName = args.option("cellSet");
|
2010-05-12 13:27:55 +00:00
|
|
|
vtkName = cellSetName;
|
|
|
|
}
|
|
|
|
else if (Pstream::parRun())
|
|
|
|
{
|
|
|
|
// Strip off leading casename, leaving just processor_DDD ending.
|
|
|
|
vtkName = runTime.caseName();
|
|
|
|
|
|
|
|
string::size_type i = vtkName.rfind("processor");
|
|
|
|
|
|
|
|
if (i != string::npos)
|
|
|
|
{
|
|
|
|
vtkName = vtkName.substr(i);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
vtkName = runTime.caseName();
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
instantList timeDirs = timeSelector::select0(runTime, args);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
# include "createNamedMesh.H"
|
|
|
|
|
|
|
|
// VTK/ directory in the case
|
|
|
|
fileName fvPath(runTime.path()/"VTK");
|
|
|
|
// Directory of mesh (region0 gets filtered out)
|
|
|
|
fileName regionPrefix = "";
|
|
|
|
|
|
|
|
if (regionName != polyMesh::defaultRegion)
|
|
|
|
{
|
|
|
|
fvPath = fvPath/regionName;
|
|
|
|
regionPrefix = regionName;
|
|
|
|
}
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
if (isDir(fvPath))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
if
|
|
|
|
(
|
2010-08-26 14:22:03 +00:00
|
|
|
args.optionFound("time")
|
|
|
|
|| args.optionFound("latestTime")
|
|
|
|
|| cellSetName.size()
|
2010-05-12 13:27:55 +00:00
|
|
|
|| regionName != polyMesh::defaultRegion
|
|
|
|
)
|
|
|
|
{
|
|
|
|
Info<< "Keeping old VTK files in " << fvPath << nl << endl;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
Info<< "Deleting old VTK files in " << fvPath << nl << endl;
|
|
|
|
|
|
|
|
rmDir(fvPath);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
mkDir(fvPath);
|
|
|
|
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
// Mesh wrapper; does subsetting and decomposition
|
2010-05-12 13:27:55 +00:00
|
|
|
vtkMesh vMesh
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
regionName,
|
|
|
|
runTime.timeName(),
|
|
|
|
runTime,
|
|
|
|
IOobject::MUST_READ
|
|
|
|
),
|
|
|
|
cellSetName
|
|
|
|
);
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
// Scan for all possible lagrangian clouds
|
|
|
|
HashSet<fileName> allCloudDirs;
|
|
|
|
forAll(timeDirs, timeI)
|
|
|
|
{
|
|
|
|
runTime.setTime(timeDirs[timeI], timeI);
|
|
|
|
fileNameList cloudDirs
|
|
|
|
(
|
|
|
|
readDir
|
|
|
|
(
|
|
|
|
runTime.timePath()/regionPrefix/cloud::prefix,
|
|
|
|
fileName::DIRECTORY
|
|
|
|
)
|
|
|
|
);
|
|
|
|
forAll(cloudDirs, i)
|
|
|
|
{
|
|
|
|
IOobjectList sprayObjs
|
|
|
|
(
|
|
|
|
mesh,
|
|
|
|
runTime.timeName(),
|
|
|
|
cloud::prefix/cloudDirs[i]
|
|
|
|
);
|
|
|
|
|
|
|
|
IOobject* positionsPtr = sprayObjs.lookup("positions");
|
|
|
|
|
|
|
|
if (positionsPtr)
|
|
|
|
{
|
|
|
|
if (allCloudDirs.insert(cloudDirs[i]))
|
|
|
|
{
|
|
|
|
Info<< "At time: " << runTime.timeName()
|
|
|
|
<< " detected cloud directory : " << cloudDirs[i]
|
|
|
|
<< endl;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
forAll(timeDirs, timeI)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
runTime.setTime(timeDirs[timeI], timeI);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
Info<< "Time: " << runTime.timeName() << endl;
|
|
|
|
|
|
|
|
word timeDesc = "";
|
|
|
|
if (useTimeName)
|
|
|
|
{
|
|
|
|
timeDesc = runTime.timeName();
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
timeDesc = name(runTime.timeIndex());
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Check for new polyMesh/ and update mesh, fvMeshSubset and cell
|
|
|
|
// decomposition.
|
|
|
|
polyMesh::readUpdateState meshState = vMesh.readUpdate();
|
|
|
|
|
|
|
|
const fvMesh& mesh = vMesh.mesh();
|
|
|
|
|
|
|
|
if
|
|
|
|
(
|
|
|
|
meshState == polyMesh::TOPO_CHANGE
|
|
|
|
|| meshState == polyMesh::TOPO_PATCH_CHANGE
|
|
|
|
)
|
|
|
|
{
|
|
|
|
Info<< " Read new mesh" << nl << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// If faceSet: write faceSet only (as polydata)
|
2010-08-26 14:22:03 +00:00
|
|
|
if (args.optionFound("faceSet"))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
// Load the faceSet
|
2010-08-26 14:22:03 +00:00
|
|
|
faceSet set(mesh, args.option("faceSet"));
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Filename as if patch with same name.
|
|
|
|
mkDir(fvPath/set.name());
|
|
|
|
|
|
|
|
fileName patchFileName
|
|
|
|
(
|
|
|
|
fvPath/set.name()/set.name()
|
|
|
|
+ "_"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk"
|
|
|
|
);
|
|
|
|
|
|
|
|
Info<< " FaceSet : " << patchFileName << endl;
|
|
|
|
|
|
|
|
writeFaceSet(binary, vMesh, set, patchFileName);
|
|
|
|
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
// If pointSet: write pointSet only (as polydata)
|
2010-08-26 14:22:03 +00:00
|
|
|
if (args.optionFound("pointSet"))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
// Load the pointSet
|
2010-08-26 14:22:03 +00:00
|
|
|
pointSet set(mesh, args.option("pointSet"));
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Filename as if patch with same name.
|
|
|
|
mkDir(fvPath/set.name());
|
|
|
|
|
|
|
|
fileName patchFileName
|
|
|
|
(
|
|
|
|
fvPath/set.name()/set.name()
|
|
|
|
+ "_"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk"
|
|
|
|
);
|
|
|
|
|
|
|
|
Info<< " pointSet : " << patchFileName << endl;
|
|
|
|
|
|
|
|
writePointSet(binary, vMesh, set, patchFileName);
|
|
|
|
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// Search for list of objects for this time
|
|
|
|
IOobjectList objects(mesh, runTime.timeName());
|
|
|
|
|
|
|
|
HashSet<word> selectedFields;
|
2010-08-26 14:22:03 +00:00
|
|
|
if (args.optionFound("fields"))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
args.optionLookup("fields")() >> selectedFields;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
// Construct the vol fields (on the original mesh if subsetted)
|
|
|
|
|
|
|
|
PtrList<volScalarField> vsf;
|
2010-09-22 18:13:13 +00:00
|
|
|
readFields(vMesh, vMesh, objects, selectedFields, vsf);
|
2010-05-12 13:27:55 +00:00
|
|
|
readFields(vMesh, vMesh, objects, selectedFields, vsf);
|
|
|
|
print(" volScalarFields :", Info, vsf);
|
|
|
|
|
|
|
|
PtrList<volVectorField> vvf;
|
|
|
|
readFields(vMesh, vMesh, objects, selectedFields, vvf);
|
|
|
|
print(" volVectorFields :", Info, vvf);
|
|
|
|
|
|
|
|
PtrList<volSphericalTensorField> vSpheretf;
|
|
|
|
readFields(vMesh, vMesh, objects, selectedFields, vSpheretf);
|
|
|
|
print(" volSphericalTensorFields :", Info, vSpheretf);
|
|
|
|
|
|
|
|
PtrList<volSymmTensorField> vSymmtf;
|
|
|
|
readFields(vMesh, vMesh, objects, selectedFields, vSymmtf);
|
|
|
|
print(" volSymmTensorFields :", Info, vSymmtf);
|
|
|
|
|
|
|
|
PtrList<volTensorField> vtf;
|
|
|
|
readFields(vMesh, vMesh, objects, selectedFields, vtf);
|
|
|
|
print(" volTensorFields :", Info, vtf);
|
|
|
|
|
|
|
|
const label nVolFields =
|
|
|
|
vsf.size()
|
|
|
|
+ vvf.size()
|
|
|
|
+ vSpheretf.size()
|
|
|
|
+ vSymmtf.size()
|
|
|
|
+ vtf.size();
|
|
|
|
|
|
|
|
|
|
|
|
// Construct pointMesh only if nessecary since constructs edge
|
|
|
|
// addressing (expensive on polyhedral meshes)
|
|
|
|
if (noPointValues)
|
|
|
|
{
|
|
|
|
Info<< " pointScalarFields : switched off"
|
|
|
|
<< " (\"-noPointValues\" option)\n";
|
|
|
|
Info<< " pointVectorFields : switched off"
|
|
|
|
<< " (\"-noPointValues\" option)\n";
|
|
|
|
}
|
|
|
|
|
|
|
|
PtrList<pointScalarField> psf;
|
|
|
|
PtrList<pointVectorField> pvf;
|
|
|
|
PtrList<pointSphericalTensorField> pSpheretf;
|
|
|
|
PtrList<pointSymmTensorField> pSymmtf;
|
|
|
|
PtrList<pointTensorField> ptf;
|
|
|
|
|
|
|
|
if (!noPointValues)
|
|
|
|
{
|
|
|
|
readFields
|
|
|
|
(
|
|
|
|
vMesh,
|
2010-09-22 18:13:13 +00:00
|
|
|
pointMesh::New(vMesh),
|
2010-05-12 13:27:55 +00:00
|
|
|
objects,
|
|
|
|
selectedFields,
|
|
|
|
psf
|
|
|
|
);
|
|
|
|
print(" pointScalarFields :", Info, psf);
|
|
|
|
|
|
|
|
readFields
|
|
|
|
(
|
|
|
|
vMesh,
|
2010-09-22 18:13:13 +00:00
|
|
|
pointMesh::New(vMesh),
|
2010-05-12 13:27:55 +00:00
|
|
|
objects,
|
|
|
|
selectedFields,
|
|
|
|
pvf
|
|
|
|
);
|
|
|
|
print(" pointVectorFields :", Info, pvf);
|
|
|
|
|
|
|
|
readFields
|
|
|
|
(
|
|
|
|
vMesh,
|
2010-09-22 18:13:13 +00:00
|
|
|
pointMesh::New(vMesh),
|
2010-05-12 13:27:55 +00:00
|
|
|
objects,
|
|
|
|
selectedFields,
|
|
|
|
pSpheretf
|
|
|
|
);
|
|
|
|
print(" pointSphericalTensorFields :", Info, pSpheretf);
|
|
|
|
|
|
|
|
readFields
|
|
|
|
(
|
|
|
|
vMesh,
|
2010-09-22 18:13:13 +00:00
|
|
|
pointMesh::New(vMesh),
|
2010-05-12 13:27:55 +00:00
|
|
|
objects,
|
|
|
|
selectedFields,
|
|
|
|
pSymmtf
|
|
|
|
);
|
|
|
|
print(" pointSymmTensorFields :", Info, pSymmtf);
|
|
|
|
|
|
|
|
readFields
|
|
|
|
(
|
|
|
|
vMesh,
|
2010-09-22 18:13:13 +00:00
|
|
|
pointMesh::New(vMesh),
|
2010-05-12 13:27:55 +00:00
|
|
|
objects,
|
|
|
|
selectedFields,
|
|
|
|
ptf
|
|
|
|
);
|
|
|
|
print(" pointTensorFields :", Info, ptf);
|
|
|
|
}
|
|
|
|
Info<< endl;
|
|
|
|
|
|
|
|
label nPointFields =
|
|
|
|
psf.size()
|
|
|
|
+ pvf.size()
|
|
|
|
+ pSpheretf.size()
|
|
|
|
+ pSymmtf.size()
|
|
|
|
+ ptf.size();
|
|
|
|
|
|
|
|
if (doWriteInternal)
|
|
|
|
{
|
|
|
|
//
|
|
|
|
// Create file and write header
|
|
|
|
//
|
|
|
|
fileName vtkFileName
|
|
|
|
(
|
|
|
|
fvPath/vtkName
|
|
|
|
+ "_"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk"
|
|
|
|
);
|
|
|
|
|
|
|
|
Info<< " Internal : " << vtkFileName << endl;
|
|
|
|
|
|
|
|
// Write mesh
|
|
|
|
internalWriter writer(vMesh, binary, vtkFileName);
|
|
|
|
|
|
|
|
// VolFields + cellID
|
|
|
|
writeFuns::writeCellDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
vMesh.nFieldCells(),
|
|
|
|
1+nVolFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write cellID field
|
|
|
|
writer.writeCellIDs();
|
|
|
|
|
|
|
|
// Write volFields
|
|
|
|
writer.write(vsf);
|
|
|
|
writer.write(vvf);
|
|
|
|
writer.write(vSpheretf);
|
|
|
|
writer.write(vSymmtf);
|
|
|
|
writer.write(vtf);
|
|
|
|
|
|
|
|
if (!noPointValues)
|
|
|
|
{
|
|
|
|
writeFuns::writePointDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
vMesh.nFieldPoints(),
|
|
|
|
nVolFields+nPointFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// pointFields
|
|
|
|
writer.write(psf);
|
|
|
|
writer.write(pvf);
|
|
|
|
writer.write(pSpheretf);
|
|
|
|
writer.write(pSymmtf);
|
|
|
|
writer.write(ptf);
|
|
|
|
|
|
|
|
// Interpolated volFields
|
2010-08-26 14:22:03 +00:00
|
|
|
volPointInterpolation pInterp(mesh);
|
2010-05-12 13:27:55 +00:00
|
|
|
writer.write(pInterp, vsf);
|
|
|
|
writer.write(pInterp, vvf);
|
|
|
|
writer.write(pInterp, vSpheretf);
|
|
|
|
writer.write(pInterp, vSymmtf);
|
|
|
|
writer.write(pInterp, vtf);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
//
|
|
|
|
// Write surface fields
|
|
|
|
//
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
if (args.optionFound("surfaceFields"))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
PtrList<surfaceScalarField> ssf;
|
|
|
|
readFields
|
|
|
|
(
|
|
|
|
vMesh,
|
|
|
|
vMesh,
|
|
|
|
objects,
|
|
|
|
selectedFields,
|
|
|
|
ssf
|
|
|
|
);
|
|
|
|
print(" surfScalarFields :", Info, ssf);
|
|
|
|
|
|
|
|
PtrList<surfaceVectorField> svf;
|
|
|
|
readFields
|
|
|
|
(
|
|
|
|
vMesh,
|
|
|
|
vMesh,
|
|
|
|
objects,
|
|
|
|
selectedFields,
|
|
|
|
svf
|
|
|
|
);
|
|
|
|
print(" surfVectorFields :", Info, svf);
|
|
|
|
|
|
|
|
if (ssf.size() + svf.size() > 0)
|
|
|
|
{
|
|
|
|
// Rework the scalar fields into vectorfields.
|
|
|
|
label sz = svf.size();
|
|
|
|
|
|
|
|
svf.setSize(sz+ssf.size());
|
|
|
|
|
|
|
|
surfaceVectorField n(mesh.Sf()/mesh.magSf());
|
|
|
|
|
|
|
|
forAll(ssf, i)
|
|
|
|
{
|
|
|
|
svf.set(sz+i, ssf[i]*n);
|
|
|
|
svf[sz+i].rename(ssf[i].name());
|
|
|
|
}
|
|
|
|
ssf.clear();
|
|
|
|
|
|
|
|
mkDir(fvPath / "surfaceFields");
|
|
|
|
|
|
|
|
fileName surfFileName
|
|
|
|
(
|
|
|
|
fvPath
|
|
|
|
/"surfaceFields"
|
|
|
|
/"surfaceFields"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ "_"
|
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk"
|
|
|
|
);
|
|
|
|
|
|
|
|
writeSurfFields
|
|
|
|
(
|
|
|
|
binary,
|
|
|
|
vMesh,
|
|
|
|
surfFileName,
|
|
|
|
svf
|
|
|
|
);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
//
|
|
|
|
// Write patches (POLYDATA file, one for each patch)
|
|
|
|
//
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
|
|
|
|
const polyBoundaryMesh& patches = mesh.boundaryMesh();
|
|
|
|
|
|
|
|
if (allPatches)
|
|
|
|
{
|
|
|
|
mkDir(fvPath/"allPatches");
|
|
|
|
|
|
|
|
fileName patchFileName;
|
|
|
|
|
|
|
|
if (vMesh.useSubMesh())
|
|
|
|
{
|
|
|
|
patchFileName =
|
|
|
|
fvPath/"allPatches"/cellSetName
|
|
|
|
+ "_"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk";
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
patchFileName =
|
|
|
|
fvPath/"allPatches"/"allPatches"
|
|
|
|
+ "_"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk";
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< " Combined patches : " << patchFileName << endl;
|
|
|
|
|
|
|
|
patchWriter writer
|
|
|
|
(
|
|
|
|
vMesh,
|
|
|
|
binary,
|
|
|
|
nearCellValue,
|
|
|
|
patchFileName,
|
|
|
|
getSelectedPatches(patches, excludePatches)
|
|
|
|
);
|
|
|
|
|
|
|
|
// VolFields + patchID
|
|
|
|
writeFuns::writeCellDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
writer.nFaces(),
|
|
|
|
1+nVolFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write patchID field
|
|
|
|
writer.writePatchIDs();
|
|
|
|
|
|
|
|
// Write volFields
|
|
|
|
writer.write(vsf);
|
|
|
|
writer.write(vvf);
|
|
|
|
writer.write(vSpheretf);
|
|
|
|
writer.write(vSymmtf);
|
|
|
|
writer.write(vtf);
|
|
|
|
|
|
|
|
if (!noPointValues)
|
|
|
|
{
|
|
|
|
writeFuns::writePointDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
writer.nPoints(),
|
2011-02-25 15:25:43 +00:00
|
|
|
nPointFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write pointFields
|
|
|
|
writer.write(psf);
|
|
|
|
writer.write(pvf);
|
|
|
|
writer.write(pSpheretf);
|
|
|
|
writer.write(pSymmtf);
|
|
|
|
writer.write(ptf);
|
|
|
|
|
|
|
|
// no interpolated volFields since I cannot be bothered to
|
|
|
|
// create the patchInterpolation for all subpatches.
|
|
|
|
}
|
|
|
|
}
|
|
|
|
else if (allWallPatches)
|
|
|
|
{
|
|
|
|
mkDir(fvPath/"allWallPatches");
|
|
|
|
|
|
|
|
fileName wallPatchFileName;
|
|
|
|
|
|
|
|
if (vMesh.useSubMesh())
|
|
|
|
{
|
|
|
|
wallPatchFileName =
|
|
|
|
fvPath/"allWallPatches"/cellSetName
|
|
|
|
+ "_"
|
|
|
|
+ timeDesc
|
|
|
|
+ ".vtk";
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
wallPatchFileName =
|
|
|
|
fvPath/"allWallPatches"/"allWallPatches"
|
|
|
|
+ "_"
|
|
|
|
+ timeDesc
|
|
|
|
+ ".vtk";
|
|
|
|
}
|
|
|
|
|
|
|
|
// Go through the boundary and add all patches that are not
|
|
|
|
// of type wall onto the exclude{atches list
|
|
|
|
forAll (patches, patchI)
|
|
|
|
{
|
|
|
|
if (!isA<wallPolyPatch>(patches[patchI]))
|
|
|
|
{
|
|
|
|
excludePatches.insert(patches[patchI].name());
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< " Combined wall patches : " << wallPatchFileName
|
|
|
|
<< endl;
|
|
|
|
|
|
|
|
patchWriter writer
|
|
|
|
(
|
|
|
|
vMesh,
|
|
|
|
binary,
|
|
|
|
nearCellValue,
|
|
|
|
wallPatchFileName,
|
|
|
|
getSelectedPatches(patches, excludePatches)
|
|
|
|
);
|
|
|
|
|
|
|
|
// VolFields + patchID
|
|
|
|
writeFuns::writeCellDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
writer.nFaces(),
|
|
|
|
1+nVolFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write patchID field
|
|
|
|
writer.writePatchIDs();
|
|
|
|
|
|
|
|
// Write volFields
|
|
|
|
writer.write(vsf);
|
|
|
|
writer.write(vvf);
|
|
|
|
writer.write(vSpheretf);
|
|
|
|
writer.write(vSymmtf);
|
|
|
|
writer.write(vtf);
|
|
|
|
|
|
|
|
if (!noPointValues)
|
|
|
|
{
|
|
|
|
writeFuns::writePointDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
writer.nPoints(),
|
2010-05-12 13:27:55 +00:00
|
|
|
nPointFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write pointFields
|
|
|
|
writer.write(psf);
|
|
|
|
writer.write(pvf);
|
|
|
|
writer.write(pSpheretf);
|
|
|
|
writer.write(pSymmtf);
|
|
|
|
writer.write(ptf);
|
|
|
|
|
|
|
|
// no interpolated volFields since I cannot be bothered to
|
|
|
|
// create the patchInterpolation for all subpatches.
|
|
|
|
}
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
forAll(patches, patchI)
|
|
|
|
{
|
|
|
|
const polyPatch& pp = patches[patchI];
|
|
|
|
|
|
|
|
if (!excludePatches.found(pp.name()))
|
|
|
|
{
|
|
|
|
mkDir(fvPath/pp.name());
|
|
|
|
|
|
|
|
fileName patchFileName;
|
|
|
|
|
|
|
|
if (vMesh.useSubMesh())
|
|
|
|
{
|
|
|
|
patchFileName =
|
|
|
|
fvPath/pp.name()/cellSetName
|
|
|
|
+ "_"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk";
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
patchFileName =
|
|
|
|
fvPath/pp.name()/pp.name()
|
|
|
|
+ "_"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk";
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< " Patch : " << patchFileName << endl;
|
|
|
|
|
|
|
|
patchWriter writer
|
|
|
|
(
|
|
|
|
vMesh,
|
|
|
|
binary,
|
|
|
|
nearCellValue,
|
|
|
|
patchFileName,
|
|
|
|
labelList(1, patchI)
|
|
|
|
);
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
if (!isA<emptyPolyPatch>(pp))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
// VolFields + patchID
|
|
|
|
writeFuns::writeCellDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
writer.nFaces(),
|
|
|
|
1+nVolFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write patchID field
|
|
|
|
writer.writePatchIDs();
|
|
|
|
|
|
|
|
// Write volFields
|
|
|
|
writer.write(vsf);
|
|
|
|
writer.write(vvf);
|
|
|
|
writer.write(vSpheretf);
|
|
|
|
writer.write(vSymmtf);
|
|
|
|
writer.write(vtf);
|
|
|
|
|
|
|
|
if (!noPointValues)
|
|
|
|
{
|
|
|
|
writeFuns::writePointDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
writer.nPoints(),
|
|
|
|
nVolFields
|
|
|
|
+ nPointFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write pointFields
|
|
|
|
writer.write(psf);
|
|
|
|
writer.write(pvf);
|
|
|
|
writer.write(pSpheretf);
|
|
|
|
writer.write(pSymmtf);
|
|
|
|
writer.write(ptf);
|
|
|
|
|
|
|
|
PrimitivePatchInterpolation<primitivePatch> pInter
|
|
|
|
(
|
|
|
|
pp
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write interpolated volFields
|
|
|
|
writer.write(pInter, vsf);
|
|
|
|
writer.write(pInter, vvf);
|
|
|
|
writer.write(pInter, vSpheretf);
|
|
|
|
writer.write(pInter, vSymmtf);
|
|
|
|
writer.write(pInter, vtf);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
//
|
|
|
|
// Write faceZones (POLYDATA file, one for each zone)
|
|
|
|
//
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
|
|
|
|
if (doFaceZones)
|
|
|
|
{
|
|
|
|
const faceZoneMesh& zones = mesh.faceZones();
|
|
|
|
|
|
|
|
forAll(zones, zoneI)
|
|
|
|
{
|
|
|
|
const faceZone& pp = zones[zoneI];
|
|
|
|
|
|
|
|
mkDir(fvPath/pp.name());
|
|
|
|
|
|
|
|
fileName patchFileName;
|
|
|
|
|
|
|
|
if (vMesh.useSubMesh())
|
|
|
|
{
|
|
|
|
patchFileName =
|
|
|
|
fvPath/pp.name()/cellSetName
|
|
|
|
+ "_"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk";
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
patchFileName =
|
|
|
|
fvPath/pp.name()/pp.name()
|
|
|
|
+ "_"
|
2010-08-26 14:22:03 +00:00
|
|
|
+ timeDesc
|
2010-05-12 13:27:55 +00:00
|
|
|
+ ".vtk";
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< " FaceZone : " << patchFileName << endl;
|
|
|
|
|
|
|
|
std::ofstream str(patchFileName.c_str());
|
|
|
|
|
|
|
|
writeFuns::writeHeader(str, binary, pp.name());
|
|
|
|
str << "DATASET POLYDATA" << std::endl;
|
|
|
|
|
|
|
|
writePatchGeom
|
|
|
|
(
|
|
|
|
binary,
|
|
|
|
pp().localFaces(),
|
|
|
|
pp().localPoints(),
|
|
|
|
str
|
|
|
|
);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
//
|
2010-09-22 18:13:13 +00:00
|
|
|
// Write finite area data
|
2010-05-12 13:27:55 +00:00
|
|
|
//
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
|
2010-09-22 18:13:13 +00:00
|
|
|
if (args.options().found("faMesh"))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2010-09-22 18:13:13 +00:00
|
|
|
mkDir(fvPath/"faMesh");
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2010-09-22 18:13:13 +00:00
|
|
|
fileName faFileName =
|
|
|
|
fvPath/"faMesh"/"faMesh"
|
|
|
|
+ "_"
|
|
|
|
+ name(timeI)
|
|
|
|
+ ".vtk";
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Create FA mesh
|
|
|
|
faMesh aMesh(vMesh.mesh());
|
|
|
|
|
|
|
|
PtrList<areaScalarField> asf;
|
|
|
|
readFieldsNoSubset(vMesh, aMesh, objects, selectedFields, asf);
|
|
|
|
print(" areaScalarFields :", Info, asf);
|
|
|
|
|
|
|
|
PtrList<areaVectorField> avf;
|
|
|
|
readFieldsNoSubset(vMesh, aMesh, objects, selectedFields, avf);
|
|
|
|
print(" areaVectorFields :", Info, avf);
|
|
|
|
|
|
|
|
const label nAreaFields =
|
|
|
|
asf.size()
|
|
|
|
+ avf.size();
|
|
|
|
|
|
|
|
faMeshWriter writer
|
|
|
|
(
|
|
|
|
aMesh,
|
|
|
|
binary,
|
|
|
|
faFileName
|
|
|
|
);
|
|
|
|
|
|
|
|
writeFuns::writeCellDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
aMesh.nFaces(),
|
|
|
|
nAreaFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write areaFields
|
|
|
|
writer.write(asf);
|
|
|
|
writer.write(avf);
|
|
|
|
|
|
|
|
if (!noPointValues)
|
|
|
|
{
|
|
|
|
writeFuns::writePointDataHeader
|
|
|
|
(
|
|
|
|
writer.os(),
|
|
|
|
aMesh.nPoints(),
|
|
|
|
nAreaFields
|
|
|
|
);
|
|
|
|
|
|
|
|
// Interpolated areaFields
|
|
|
|
PrimitivePatchInterpolation<indirectPrimitivePatch> pInterp
|
|
|
|
(
|
|
|
|
aMesh.patch()
|
|
|
|
);
|
|
|
|
writer.write(pInterp, asf);
|
|
|
|
writer.write(pInterp, avf);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
//
|
|
|
|
// Write lagrangian data
|
|
|
|
//
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
|
2010-09-22 18:13:13 +00:00
|
|
|
forAllConstIter(HashSet<fileName>, allCloudDirs, iter)
|
|
|
|
{
|
|
|
|
const fileName& cloudName = iter.key();
|
|
|
|
|
|
|
|
// Always create the cloud directory.
|
|
|
|
mkDir(fvPath/cloud::prefix/cloudName);
|
|
|
|
|
|
|
|
fileName lagrFileName
|
|
|
|
(
|
|
|
|
fvPath/cloud::prefix/cloudName/cloudName
|
|
|
|
+ "_" + timeDesc + ".vtk"
|
|
|
|
);
|
|
|
|
|
|
|
|
Info<< " Lagrangian: " << lagrFileName << endl;
|
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
IOobjectList sprayObjs
|
|
|
|
(
|
|
|
|
mesh,
|
|
|
|
runTime.timeName(),
|
2010-08-26 14:22:03 +00:00
|
|
|
cloud::prefix/cloudName
|
2010-05-12 13:27:55 +00:00
|
|
|
);
|
|
|
|
|
|
|
|
IOobject* positionsPtr = sprayObjs.lookup("positions");
|
|
|
|
|
|
|
|
if (positionsPtr)
|
|
|
|
{
|
|
|
|
wordList labelNames(sprayObjs.names(labelIOField::typeName));
|
|
|
|
Info<< " labels :";
|
|
|
|
print(Info, labelNames);
|
|
|
|
|
|
|
|
wordList scalarNames(sprayObjs.names(scalarIOField::typeName));
|
|
|
|
Info<< " scalars :";
|
|
|
|
print(Info, scalarNames);
|
|
|
|
|
|
|
|
wordList vectorNames(sprayObjs.names(vectorIOField::typeName));
|
|
|
|
Info<< " vectors :";
|
|
|
|
print(Info, vectorNames);
|
|
|
|
|
|
|
|
wordList sphereNames
|
|
|
|
(
|
|
|
|
sprayObjs.names
|
|
|
|
(
|
|
|
|
sphericalTensorIOField::typeName
|
|
|
|
)
|
|
|
|
);
|
|
|
|
Info<< " spherical tensors :";
|
|
|
|
print(Info, sphereNames);
|
|
|
|
|
|
|
|
wordList symmNames
|
|
|
|
(
|
|
|
|
sprayObjs.names
|
|
|
|
(
|
|
|
|
symmTensorIOField::typeName
|
|
|
|
)
|
|
|
|
);
|
|
|
|
Info<< " symm tensors :";
|
|
|
|
print(Info, symmNames);
|
|
|
|
|
|
|
|
wordList tensorNames(sprayObjs.names(tensorIOField::typeName));
|
|
|
|
Info<< " tensors :";
|
|
|
|
print(Info, tensorNames);
|
|
|
|
|
|
|
|
lagrangianWriter writer
|
|
|
|
(
|
|
|
|
vMesh,
|
|
|
|
binary,
|
|
|
|
lagrFileName,
|
2010-08-26 14:22:03 +00:00
|
|
|
cloudName,
|
|
|
|
false
|
2010-05-12 13:27:55 +00:00
|
|
|
);
|
|
|
|
|
|
|
|
// Write number of fields
|
|
|
|
writer.writeParcelHeader
|
|
|
|
(
|
|
|
|
labelNames.size()
|
2010-08-26 14:22:03 +00:00
|
|
|
+ scalarNames.size()
|
|
|
|
+ vectorNames.size()
|
|
|
|
+ sphereNames.size()
|
|
|
|
+ symmNames.size()
|
|
|
|
+ tensorNames.size()
|
2010-05-12 13:27:55 +00:00
|
|
|
);
|
|
|
|
|
|
|
|
// Fields
|
|
|
|
writer.writeIOField<label>(labelNames);
|
|
|
|
writer.writeIOField<scalar>(scalarNames);
|
|
|
|
writer.writeIOField<vector>(vectorNames);
|
2010-08-26 14:22:03 +00:00
|
|
|
writer.writeIOField<sphericalTensor>(sphereNames);
|
|
|
|
writer.writeIOField<symmTensor>(symmNames);
|
|
|
|
writer.writeIOField<tensor>(tensorNames);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
lagrangianWriter writer
|
|
|
|
(
|
|
|
|
vMesh,
|
|
|
|
binary,
|
|
|
|
lagrFileName,
|
|
|
|
cloudName,
|
|
|
|
true
|
|
|
|
);
|
|
|
|
|
|
|
|
// Write number of fields
|
|
|
|
writer.writeParcelHeader(0);
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
//
|
|
|
|
// Link parallel outputs back to undecomposed case for ease of loading
|
|
|
|
//
|
|
|
|
//---------------------------------------------------------------------
|
|
|
|
|
|
|
|
if (Pstream::parRun() && doLinks)
|
|
|
|
{
|
|
|
|
mkDir(runTime.path()/".."/"VTK");
|
|
|
|
chDir(runTime.path()/".."/"VTK");
|
|
|
|
|
|
|
|
Info<< "Linking all processor files to " << runTime.path()/".."/"VTK"
|
|
|
|
<< endl;
|
|
|
|
|
|
|
|
// Get list of vtk files
|
|
|
|
fileName procVTK
|
|
|
|
(
|
|
|
|
fileName("..")
|
|
|
|
/"processor" + name(Pstream::myProcNo())
|
|
|
|
/"VTK"
|
|
|
|
);
|
|
|
|
|
|
|
|
fileNameList dirs(readDir(procVTK, fileName::DIRECTORY));
|
|
|
|
label sz = dirs.size();
|
|
|
|
dirs.setSize(sz+1);
|
|
|
|
dirs[sz] = ".";
|
|
|
|
|
|
|
|
forAll(dirs, i)
|
|
|
|
{
|
|
|
|
fileNameList subFiles(readDir(procVTK/dirs[i], fileName::FILE));
|
|
|
|
|
|
|
|
forAll(subFiles, j)
|
|
|
|
{
|
|
|
|
fileName procFile(procVTK/dirs[i]/subFiles[j]);
|
|
|
|
|
|
|
|
if (exists(procFile))
|
|
|
|
{
|
|
|
|
string cmd
|
|
|
|
(
|
|
|
|
"ln -s "
|
|
|
|
+ procFile
|
|
|
|
+ " "
|
|
|
|
+ "processor"
|
|
|
|
+ name(Pstream::myProcNo())
|
|
|
|
+ "_"
|
|
|
|
+ procFile.name()
|
|
|
|
);
|
|
|
|
system(cmd.c_str());
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
|
|
|
|
return 0;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|