createTurbulenceFields: Specification of turbulence fields now command-line option

Usage: createTurbulenceFields [OPTIONS]
options:
  -case <dir>       specify alternate case directory, default is the cwd
  -constant         include the 'constant/' dir in the times list
  -fields <wordReList>
                    specify which turbulence fields (k, epsilon, omega, R) to
                    write - eg '(k omega)' or '(R)' or '(.*)'.
  -latestTime       select the latest time
  -newTimes         select the new times
  -noFunctionObjects
                    do not execute functionObjects
  -noZero           exclude the '0/' dir from the times list, has precedence
                    over the -withZero option
  -parallel         run in parallel
  -roots <(dir1 .. dirN)>
                    slave root directories for distributed running
  -time <ranges>    comma-separated time ranges - eg, ':10,20,40:70,1000:'
  -srcDoc           display source code in browser
  -doc              display application documentation in browser
  -help             print the usage

Resolves feature request http://www.openfoam.org/mantisbt/view.php?id=1912
This commit is contained in:
Henry Weller
2015-11-13 10:55:34 +00:00
parent 778c0cc8ed
commit c18fd1d300

View File

@ -44,11 +44,26 @@ int main(int argc, char *argv[])
{ {
timeSelector::addOptions(); timeSelector::addOptions();
argList::addOption
(
"fields",
"wordReList",
"specify which turbulence fields (k, epsilon, omega, R) to write"
" - eg '(k omega)' or '(R)' or '(.*)'."
);
#include "setRootCase.H" #include "setRootCase.H"
#include "createTime.H" #include "createTime.H"
instantList timeDirs = timeSelector::select0(runTime, args); instantList timeDirs = timeSelector::select0(runTime, args);
const bool selectedFields = args.optionFound("fields");
wordReList fieldPatterns;
if (selectedFields)
{
fieldPatterns = wordReList(args.optionLookup("fields")());
}
#include "createMesh.H" #include "createMesh.H"
forAll(timeDirs, timeI) forAll(timeDirs, timeI)
@ -59,69 +74,77 @@ int main(int argc, char *argv[])
#include "createFields.H" #include "createFields.H"
Info<< "\nRetrieving field k from turbulence model" << endl; if (findStrings(fieldPatterns, "k"))
const volScalarField k(RASModel->k());
Info<< "\nRetrieving field epsilon from turbulence model" << endl;
const volScalarField epsilon(RASModel->epsilon());
Info<< "\nRetrieving field R from turbulence model" << endl;
const volSymmTensorField R(RASModel->R());
// Check availability of tubulence fields
if (!IOobject("k", runTime.timeName(), mesh).headerOk())
{ {
Info<< "\nWriting turbulence field k" << endl; if (!IOobject("k", runTime.timeName(), mesh).headerOk())
k.write(); {
} Info<< " Writing turbulence field k" << endl;
else const volScalarField k(RASModel->k());
{ k.write();
Info<< "\nTurbulence k field already exists" << endl; }
else
{
Info<< " Turbulence k field already exists" << endl;
}
} }
if (!IOobject("epsilon", runTime.timeName(), mesh).headerOk()) if (findStrings(fieldPatterns, "epsilon"))
{ {
Info<< "\nWriting turbulence field epsilon" << endl; if (!IOobject("epsilon", runTime.timeName(), mesh).headerOk())
epsilon.write(); {
} Info<< " Writing turbulence field epsilon" << endl;
else const volScalarField epsilon(RASModel->epsilon());
{ epsilon.write();
Info<< "\nTurbulence epsilon field already exists" << endl; }
else
{
Info<< " Turbulence epsilon field already exists" << endl;
}
} }
if (!IOobject("R", runTime.timeName(), mesh).headerOk()) if (findStrings(fieldPatterns, "R"))
{ {
Info<< "\nWriting turbulence field R" << endl; if (!IOobject("R", runTime.timeName(), mesh).headerOk())
R.write(); {
} Info<< " Writing turbulence field R" << endl;
else const volSymmTensorField R(RASModel->R());
{ R.write();
Info<< "\nTurbulence R field already exists" << endl; }
else
{
Info<< " Turbulence R field already exists" << endl;
}
} }
if (!IOobject("omega", runTime.timeName(), mesh).headerOk()) if (findStrings(fieldPatterns, "omega"))
{ {
const scalar Cmu = 0.09; if (!IOobject("omega", runTime.timeName(), mesh).headerOk())
{
const scalar Cmu = 0.09;
Info<< "creating omega" << endl; // Assume k and epsilon are available
volScalarField omega const volScalarField k(RASModel->k());
( const volScalarField epsilon(RASModel->epsilon());
IOobject
volScalarField omega
( (
"omega", IOobject
runTime.timeName(), (
mesh "omega",
), runTime.timeName(),
epsilon/(Cmu*k), mesh
epsilon.boundaryField().types() ),
); epsilon/(Cmu*k),
Info<< "\nWriting turbulence field omega" << endl; epsilon.boundaryField().types()
omega.write(); );
}
else Info<< " Writing turbulence field omega" << endl;
{ omega.write();
Info<< "\nTurbulence omega field already exists" << endl; }
else
{
Info<< " Turbulence omega field already exists" << endl;
}
} }
} }