mirror of
https://github.com/OpenFOAM/OpenFOAM-6.git
synced 2025-12-08 06:57:46 +00:00
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:
@ -44,11 +44,26 @@ int main(int argc, char *argv[])
|
||||
{
|
||||
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 "createTime.H"
|
||||
|
||||
instantList timeDirs = timeSelector::select0(runTime, args);
|
||||
|
||||
const bool selectedFields = args.optionFound("fields");
|
||||
wordReList fieldPatterns;
|
||||
if (selectedFields)
|
||||
{
|
||||
fieldPatterns = wordReList(args.optionLookup("fields")());
|
||||
}
|
||||
|
||||
#include "createMesh.H"
|
||||
|
||||
forAll(timeDirs, timeI)
|
||||
@ -59,69 +74,77 @@ int main(int argc, char *argv[])
|
||||
|
||||
#include "createFields.H"
|
||||
|
||||
Info<< "\nRetrieving field k from turbulence model" << endl;
|
||||
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())
|
||||
if (findStrings(fieldPatterns, "k"))
|
||||
{
|
||||
Info<< "\nWriting turbulence field k" << endl;
|
||||
k.write();
|
||||
}
|
||||
else
|
||||
{
|
||||
Info<< "\nTurbulence k field already exists" << endl;
|
||||
if (!IOobject("k", runTime.timeName(), mesh).headerOk())
|
||||
{
|
||||
Info<< " Writing turbulence field k" << endl;
|
||||
const volScalarField k(RASModel->k());
|
||||
k.write();
|
||||
}
|
||||
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;
|
||||
epsilon.write();
|
||||
}
|
||||
else
|
||||
{
|
||||
Info<< "\nTurbulence epsilon field already exists" << endl;
|
||||
if (!IOobject("epsilon", runTime.timeName(), mesh).headerOk())
|
||||
{
|
||||
Info<< " Writing turbulence field epsilon" << endl;
|
||||
const volScalarField epsilon(RASModel->epsilon());
|
||||
epsilon.write();
|
||||
}
|
||||
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;
|
||||
R.write();
|
||||
}
|
||||
else
|
||||
{
|
||||
Info<< "\nTurbulence R field already exists" << endl;
|
||||
if (!IOobject("R", runTime.timeName(), mesh).headerOk())
|
||||
{
|
||||
Info<< " Writing turbulence field R" << endl;
|
||||
const volSymmTensorField R(RASModel->R());
|
||||
R.write();
|
||||
}
|
||||
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;
|
||||
volScalarField omega
|
||||
(
|
||||
IOobject
|
||||
// Assume k and epsilon are available
|
||||
const volScalarField k(RASModel->k());
|
||||
const volScalarField epsilon(RASModel->epsilon());
|
||||
|
||||
volScalarField omega
|
||||
(
|
||||
"omega",
|
||||
runTime.timeName(),
|
||||
mesh
|
||||
),
|
||||
epsilon/(Cmu*k),
|
||||
epsilon.boundaryField().types()
|
||||
);
|
||||
Info<< "\nWriting turbulence field omega" << endl;
|
||||
omega.write();
|
||||
}
|
||||
else
|
||||
{
|
||||
Info<< "\nTurbulence omega field already exists" << endl;
|
||||
IOobject
|
||||
(
|
||||
"omega",
|
||||
runTime.timeName(),
|
||||
mesh
|
||||
),
|
||||
epsilon/(Cmu*k),
|
||||
epsilon.boundaryField().types()
|
||||
);
|
||||
|
||||
Info<< " Writing turbulence field omega" << endl;
|
||||
omega.write();
|
||||
}
|
||||
else
|
||||
{
|
||||
Info<< " Turbulence omega field already exists" << endl;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
Reference in New Issue
Block a user