Go to the documentation of this file.
93 if (nameList[i] ==
name)
104 int main(
int argc,
char *argv[])
112 "write in ASCII format instead of 'C Binary'"
117 "write values in nodes"
122 "suppress writing any patches"
128 "specify particular patches to write - eg '(outlet \"inlet.*\")'. "
129 "An empty list suppresses writing the internalMesh."
135 "specify faceZones to write - eg '( slice \"mfp-.*\" )'."
141 "specify fields to export (all by default) - eg '( \"U.*\" )'."
147 "specify cellZone to write"
170 const label nVolFieldTypes = 10;
171 const word volFieldTypes[] =
173 volScalarField::typeName,
174 volVectorField::typeName,
175 volSphericalTensorField::typeName,
176 volSymmTensorField::typeName,
177 volTensorField::typeName,
179 volScalarField::DimensionedInternalField::typeName,
180 volVectorField::DimensionedInternalField::typeName,
181 volSphericalTensorField::DimensionedInternalField::typeName,
182 volSymmTensorField::DimensionedInternalField::typeName,
183 volTensorField::DimensionedInternalField::typeName
192 if (
isDir(ensightDir))
205 OFstream *ensightCaseFilePtr = NULL;
208 fileName caseFileName = prepend +
"case";
209 Info<<
nl <<
"write case: " << caseFileName.c_str() <<
endl;
214 ensightDir/caseFileName,
220 <<
"type: ensight gold" <<
nl <<
nl;
223 OFstream& ensightCaseFile = *ensightCaseFilePtr;
252 Info<<
"Converting cellZone " << cellZoneName
253 <<
" only (puts outside faces into patch "
258 meshSubsetter.setLargeCellSubset(c0, 0);
264 meshSubsetter.hasSubMesh()
265 ? meshSubsetter.subMesh()
266 : meshSubsetter.baseMesh()
277 runTime.setTime(Times.last(), Times.
size()-1);
285 Info<<
"Detected a moving mesh (multiple polyMesh/points files)."
286 <<
" Writing meshes for every timestep." <<
endl;
293 word geomFileName = prepend +
"0000";
298 geomFileName = prepend +
"****";
304 << (geomFileName +
".mesh").c_str() <<
nl;
311 runTime.setTime(Times[timeI], timeI);
328 IOobject* positionsPtr = cloudObjs.lookup(
word(
"positions"));
332 allCloudNames.
insert(cloudDirs[cloudI]);
358 allCloudFields.
find(cloudIter.key());
364 runTime.setTime(Times[timeI], timeI);
377 if (obj.
name() !=
"positions")
390 label nTimeSteps = 0;
399 Info<<
"Translating time = " << runTime.timeName() <<
nl;
402 if (
timeIndex != 0 && meshSubsetter.hasSubMesh())
404 Info<<
"Converting cellZone " << cellZoneName
405 <<
" only (puts outside faces into patch "
410 meshSubsetter.setLargeCellSubset(c0, 0);
437 ensightCaseFile<<
nl <<
"VARIABLE" <<
nl;
444 for (
label i=0; i<nVolFieldTypes; i++)
477 if (volFieldTypes[i] == volScalarField::typeName)
492 else if (volFieldTypes[i] == volVectorField::typeName)
507 else if (volFieldTypes[i] == volSphericalTensorField::typeName)
510 ensightField<sphericalTensor>
522 else if (volFieldTypes[i] == volSymmTensorField::typeName)
525 ensightField<symmTensor>
537 else if (volFieldTypes[i] == volTensorField::typeName)
556 == volScalarField::DimensionedInternalField::typeName
566 volField<scalar>(meshSubsetter, df),
579 == volVectorField::DimensionedInternalField::typeName
589 volField<vector>(meshSubsetter, df),
602 == volSphericalTensorField::DimensionedInternalField::typeName
610 ensightField<sphericalTensor>
612 volField<sphericalTensor>(meshSubsetter, df),
625 == volSymmTensorField::DimensionedInternalField::typeName
633 ensightField<symmTensor>
635 volField<symmTensor>(meshSubsetter, df),
648 == volTensorField::DimensionedInternalField::typeName
658 volField<tensor>(meshSubsetter, df),
685 bool cloudExists = inFileNameList(currentCloudDirs,
cloudName);
697 const word& fieldName = fieldIter.key();
698 const word& fieldType = fieldIter();
712 ensightCloudField<scalar>
725 ensightCloudField<vector>
738 Info<<
"Unable to convert field type " << fieldType
739 <<
" for field " << fieldName <<
endl;
749 delete ensightCaseFilePtr;
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A class for handling words, derived from string.
A class for handling file names.
static void addOption(const word &opt, const string ¶m="", const string &usage="")
Add to an option to validOptions with usage information.
static word defaultRegion
Return the default region name.
#define forAll(list, i)
Loop across all elements in list.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
static const char *const typeName
List< wordRe > wordReList
A List of wordRe (word or regular expression)
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
const cellZoneMesh & cellZones() const
Return cell zone mesh.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
A HashTable with keys but without contents.
const fileName & rootPath() const
Return root path.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
bool headerOk()
Read and check header info.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & headerClassName() const
Return name of the class name read from header.
Istream and Ostream manipulators taking arguments.
const word & name() const
Return name.
int main(int argc, char *argv[])
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
A collection of cell labels.
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
List of IOobjects with searching and retrieving facilities.
static bool master(const label communicator=0)
Am I the master process.
readUpdateState
Enumeration defining the state of the mesh after a read update.
An STL-conforming hash table.
static const word prefix
The prefix to local: lagrangian.
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
wordList fieldNames(const IOobjectList &objects, const bool syncPar)
Get sorted names of fields of type. If syncPar and running in parallel.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
bool insert(const Key &key)
Insert a new entry.
const word cloudName(propsDict.lookup("cloudName"))
bool optionFound(const word &opt) const
Return true if the named option is found.
void ensightParticlePositions(const polyMesh &mesh, const fileName &dataDir, const fileName &subDir, const word &cloudName, IOstream::streamFormat format)
const Time & time() const
Return the top-level database.
IOobject fieldObject(fieldNames[var2field[nVar]], runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE)
void size(const label)
Override size to be inconsistent with allocated storage.
fileNameList readDir(const fileName &, const fileName::Type=fileName::FILE, const bool filtergz=true)
Read a directory and return the entries as a string list.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Generic GeometricField class.
Foam::argList args(argc, argv)
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > volField(const Foam::fvMeshSubset &, const Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > &vf)
Wrapper to get hold of the field or the subsetted field.
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
bool rmDir(const fileName &)
Remove a dirctory and its contents.
const fileName & globalCaseName() const
Return case name.
word name(const complex &)
Return a string representation of a complex.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.