1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM 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.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
foamToEnsight
Description
Translates OpenFOAM data to EnSight format.
An Ensight part is created for the internalMesh and for each patch.
Usage
\b foamToEnsight [OPTION]
Options:
- \par -ascii
Write Ensight data in ASCII format instead of "C Binary"
- \par -patches patchList
Specify particular patches to write.
Specifying an empty list suppresses writing the internalMesh.
- \par -noPatches
Suppress writing any patches.
- \par -faceZones zoneList
Specify faceZones to write, with wildcards
- \par -cellZone zoneName
Specify single cellZone to write (not lagrangian)
Note
Parallel support for cloud data is not supported
- writes to \a EnSight directory to avoid collisions with foamToEnsightParts
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "timeSelector.H"
#include "IOobjectList.H"
#include "IOmanip.H"
#include "OFstream.H"
#include "volFields.H"
#include "labelIOField.H"
#include "scalarIOField.H"
#include "tensorIOField.H"
#include "ensightMesh.H"
#include "ensightField.H"
#include "ensightParticlePositions.H"
#include "ensightCloudField.H"
#include "fvc.H"
#include "cellSet.H"
#include "fvMeshSubset.H"
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
bool inFileNameList
(
const fileNameList& nameList,
const word& name
)
{
forAll(nameList, i)
{
if (nameList[i] == name)
{
return true;
}
}
return false;
}
int main(int argc, char *argv[])
{
timeSelector::addOptions();
#include "addRegionOption.H"
argList::addBoolOption
(
"ascii",
"write in ASCII format instead of 'C Binary'"
);
argList::addBoolOption
(
"nodeValues",
"write values in nodes"
);
argList::addBoolOption
(
"noPatches",
"suppress writing any patches"
);
argList::addOption
(
"patches",
"wordReList",
"specify particular patches to write - eg '(outlet \"inlet.*\")'. "
"An empty list suppresses writing the internalMesh."
);
argList::addOption
(
"faceZones",
"wordReList",
"specify faceZones to write - eg '( slice \"mfp-.*\" )'."
);
argList::addOption
(
"fields",
"wordReList",
"specify fields to export (all by default) - eg '( \"U.*\" )'."
);
argList::addOption
(
"cellZone",
"word",
"specify cellZone to write"
);
#include "setRootCase.H"
// Check options
const bool binary = !args.optionFound("ascii");
const bool nodeValues = args.optionFound("nodeValues");
#include "createTime.H"
instantList Times = timeSelector::select0(runTime, args);
#include "createNamedMesh.H"
// Mesh instance (region0 gets filtered out)
fileName regionPrefix = "";
if (regionName != polyMesh::defaultRegion)
{
regionPrefix = regionName;
}
const label nVolFieldTypes = 5;
const word volFieldTypes[] =
{
volScalarField::typeName,
volVectorField::typeName,
volSphericalTensorField::typeName,
volSymmTensorField::typeName,
volTensorField::typeName
};
// Path to EnSight directory at case level only
// - For parallel cases, data only written from master
fileName ensightDir = args.rootPath()/args.globalCaseName()/"EnSight";
if (Pstream::master())
{
if (isDir(ensightDir))
{
rmDir(ensightDir);
}
mkDir(ensightDir);
}
// Start of case file header output
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
const word prepend = args.globalCaseName() + '.';
OFstream *ensightCaseFilePtr = NULL;
if (Pstream::master())
{
fileName caseFileName = prepend + "case";
Info<< nl << "write case: " << caseFileName.c_str() << endl;
// the case file is always ASCII
ensightCaseFilePtr = new OFstream
(
ensightDir/caseFileName,
IOstream::ASCII
);
*ensightCaseFilePtr
<< "FORMAT" << nl
<< "type: ensight gold" << nl << nl;
}
OFstream& ensightCaseFile = *ensightCaseFilePtr;
// Construct the EnSight mesh
const bool selectedPatches = args.optionFound("patches");
wordReList patchPatterns;
if (selectedPatches)
{
patchPatterns = wordReList(args.optionLookup("patches")());
}
const bool selectedZones = args.optionFound("faceZones");
wordReList zonePatterns;
if (selectedZones)
{
zonePatterns = wordReList(args.optionLookup("faceZones")());
}
const bool selectedFields = args.optionFound("fields");
wordReList fieldPatterns;
if (selectedFields)
{
fieldPatterns = wordReList(args.optionLookup("fields")());
}
word cellZoneName;
const bool doCellZone = args.optionReadIfPresent("cellZone", cellZoneName);
fvMeshSubset meshSubsetter(mesh);
if (doCellZone)
{
Info<< "Converting cellZone " << cellZoneName
<< " only (puts outside faces into patch "
<< mesh.boundaryMesh()[0].name()
<< ")" << endl;
const cellZone& cz = mesh.cellZones()[cellZoneName];
cellSet c0(mesh, "c0", labelHashSet(cz));
meshSubsetter.setLargeCellSubset(c0, 0);
}
ensightMesh eMesh
(
(
meshSubsetter.hasSubMesh()
? meshSubsetter.subMesh()
: meshSubsetter.baseMesh()
),
args.optionFound("noPatches"),
selectedPatches,
patchPatterns,
selectedZones,
zonePatterns,
binary
);
// Set Time to the last time before looking for the lagrangian objects
runTime.setTime(Times.last(), Times.size()-1);
IOobjectList objects(mesh, runTime.timeName());
#include "checkMeshMoving.H"
if (meshMoving)
{
Info<< "Detected a moving mesh (multiple polyMesh/points files)."
<< " Writing meshes for every timestep." << endl;
}
wordHashSet allCloudNames;
if (Pstream::master())
{
word geomFileName = prepend + "0000";
// test pre check variable if there is a moving mesh
if (meshMoving)
{
geomFileName = prepend + "****";
}
ensightCaseFile
<< "GEOMETRY" << nl
<< "model: 1 "
<< (geomFileName + ".mesh").c_str() << nl;
}
// Identify if lagrangian data exists at each time, and add clouds
// to the 'allCloudNames' hash set
forAll(Times, timeI)
{
runTime.setTime(Times[timeI], timeI);
fileNameList cloudDirs = readDir
(
runTime.timePath()/regionPrefix/cloud::prefix,
fileName::DIRECTORY
);
forAll(cloudDirs, cloudI)
{
IOobjectList cloudObjs
(
mesh,
runTime.timeName(),
cloud::prefix/cloudDirs[cloudI]
);
IOobject* positionsPtr = cloudObjs.lookup(word("positions"));
if (positionsPtr)
{
allCloudNames.insert(cloudDirs[cloudI]);
}
}
}
HashTable<HashTable<word>> allCloudFields;
forAllConstIter(wordHashSet, allCloudNames, cloudIter)
{
// Add the name of the cloud(s) to the case file header
if (Pstream::master())
{
ensightCaseFile
<< (
"measured: 1 "
+ prepend
+ "****."
+ cloudIter.key()
).c_str()
<< nl;
}
// Create a new hash table for each cloud
allCloudFields.insert(cloudIter.key(), HashTable<word>());
// Identify the new cloud in the hash table
HashTable<HashTable<word>>::iterator newCloudIter =
allCloudFields.find(cloudIter.key());
// Loop over all times to build list of fields and field types
// for each cloud
forAll(Times, timeI)
{
runTime.setTime(Times[timeI], timeI);
IOobjectList cloudObjs
(
mesh,
runTime.timeName(),
cloud::prefix/cloudIter.key()
);
forAllConstIter(IOobjectList, cloudObjs, fieldIter)
{
const IOobject obj = *fieldIter();
if (obj.name() != "positions")
{
// Add field and field type
newCloudIter().insert
(
obj.name(),
obj.headerClassName()
);
}
}
}
}
label nTimeSteps = 0;
forAll(Times, timeIndex)
{
nTimeSteps++;
runTime.setTime(Times[timeIndex], timeIndex);
word timeName = itoa(timeIndex);
word timeFile = prepend + timeName;
Info<< "Translating time = " << runTime.timeName() << nl;
polyMesh::readUpdateState meshState = mesh.readUpdate();
if (timeIndex != 0 && meshSubsetter.hasSubMesh())
{
Info<< "Converting cellZone " << cellZoneName
<< " only (puts outside faces into patch "
<< mesh.boundaryMesh()[0].name()
<< ")" << endl;
const cellZone& cz = mesh.cellZones()[cellZoneName];
cellSet c0(mesh, "c0", labelHashSet(cz));
meshSubsetter.setLargeCellSubset(c0, 0);
}
if (meshState != polyMesh::UNCHANGED)
{
eMesh.correct();
}
if (timeIndex == 0 || meshMoving)
{
eMesh.write
(
ensightDir,
prepend,
timeIndex,
meshMoving,
ensightCaseFile
);
}
// Start of field data output
// ~~~~~~~~~~~~~~~~~~~~~~~~~~
if (timeIndex == 0 && Pstream::master())
{
ensightCaseFile<< nl << "VARIABLE" << nl;
}
// Cell field data output
// ~~~~~~~~~~~~~~~~~~~~~~
for (label i=0; i<nVolFieldTypes; i++)
{
wordList fieldNames = objects.names(volFieldTypes[i]);
forAll(fieldNames, j)
{
const word& fieldName = fieldNames[j];
// Check if the field has to be exported
if (selectedFields)
{
if (!findStrings(fieldPatterns, fieldName))
{
continue;
}
}
#include "checkData.H"
if (!variableGood)
{
continue;
}
IOobject fieldObject
(
fieldName,
mesh.time().timeName(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
);
if (volFieldTypes[i] == volScalarField::typeName)
{
volScalarField vf(fieldObject, mesh);
ensightField<scalar>
(
volField(meshSubsetter, vf),
eMesh,
ensightDir,
prepend,
timeIndex,
binary,
nodeValues,
ensightCaseFile
);
}
else if (volFieldTypes[i] == volVectorField::typeName)
{
volVectorField vf(fieldObject, mesh);
ensightField<vector>
(
volField(meshSubsetter, vf),
eMesh,
ensightDir,
prepend,
timeIndex,
binary,
nodeValues,
ensightCaseFile
);
}
else if (volFieldTypes[i] == volSphericalTensorField::typeName)
{
volSphericalTensorField vf(fieldObject, mesh);
ensightField<sphericalTensor>
(
volField(meshSubsetter, vf),
eMesh,
ensightDir,
prepend,
timeIndex,
binary,
nodeValues,
ensightCaseFile
);
}
else if (volFieldTypes[i] == volSymmTensorField::typeName)
{
volSymmTensorField vf(fieldObject, mesh);
ensightField<symmTensor>
(
volField(meshSubsetter, vf),
eMesh,
ensightDir,
prepend,
timeIndex,
binary,
nodeValues,
ensightCaseFile
);
}
else if (volFieldTypes[i] == volTensorField::typeName)
{
volTensorField vf(fieldObject, mesh);
ensightField<tensor>
(
volField(meshSubsetter, vf),
eMesh,
ensightDir,
prepend,
timeIndex,
binary,
nodeValues,
ensightCaseFile
);
}
}
}
// Cloud field data output
// ~~~~~~~~~~~~~~~~~~~~~~~
forAllConstIter(HashTable<HashTable<word>>, allCloudFields, cloudIter)
{
const word& cloudName = cloudIter.key();
fileNameList currentCloudDirs = readDir
(
runTime.timePath()/regionPrefix/cloud::prefix,
fileName::DIRECTORY
);
bool cloudExists = inFileNameList(currentCloudDirs, cloudName);
ensightParticlePositions
(
mesh,
ensightDir,
timeFile,
cloudName,
cloudExists
);
forAllConstIter(HashTable<word>, cloudIter(), fieldIter)
{
const word& fieldName = fieldIter.key();
const word& fieldType = fieldIter();
IOobject fieldObject
(
fieldName,
mesh.time().timeName(),
cloud::prefix/cloudName,
mesh,
IOobject::MUST_READ
);
bool fieldExists = fieldObject.headerOk();
if (fieldType == scalarIOField::typeName)
{
ensightCloudField<scalar>
(
fieldObject,
ensightDir,
prepend,
timeIndex,
cloudName,
ensightCaseFile,
fieldExists
);
}
else if (fieldType == vectorIOField::typeName)
{
ensightCloudField<vector>
(
fieldObject,
ensightDir,
prepend,
timeIndex,
cloudName,
ensightCaseFile,
fieldExists
);
}
else
{
Info<< "Unable to convert field type " << fieldType
<< " for field " << fieldName << endl;
}
}
}
}
#include "ensightCaseTail.H"
if (Pstream::master())
{
delete ensightCaseFilePtr;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //
|