Provides better context for the available boundary conditions, fvOptions,
functionObjects etc. and thus returns only those available to and compatible
with the particular application.
e.g.
pimpleFoam -help
Usage: pimpleFoam [OPTIONS]
options:
-case <dir> specify alternate case directory, default is the cwd
-listFunctionObjects
List functionObjects
-listFvOptions List fvOptions
-listRegisteredSwitches
List switches registered for run-time modification
-listScalarBCs List scalar field boundary conditions (fvPatchField<scalar>)
-listSwitches List switches declared in libraries but not set in
etc/controlDict
-listTurbulenceModels
List turbulenceModels
-listUnsetSwitches
List switches declared in libraries but not set in
etc/controlDict
-listVectorBCs List vector field boundary conditions (fvPatchField<vector>)
-noFunctionObjects
do not execute functionObjects
-parallel run in parallel
-postProcess Execute functionObjects only
-roots <(dir1 .. dirN)>
slave root directories for distributed running
-srcDoc display source code in browser
-doc display application documentation in browser
-help print the usage
pimpleFoam listTurbulenceModels
pimpleFoam -listTurbulenceModels
/*---------------------------------------------------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
Build : dev-39c46019e44f
Exec : pimpleFoam -listTurbulenceModels
Date : Jun 10 2017
Time : 21:37:49
Host : "dm"
PID : 675
Case : /home/dm2/henry/OpenFOAM/OpenFOAM-dev
nProcs : 1
sigFpe : Enabling floating point exception trapping (FOAM_SIGFPE).
SetNaN : Initialising allocated memory to NaN (FOAM_SETNAN).
fileModificationChecking : Monitoring run-time modified files using timeStampMaster (fileModificationSkew 10)
allowSystemOperations : Allowing user-supplied system call operations
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Turbulence models
3
(
LES
RAS
laminar
)
RAS models
18
(
LRR
LamBremhorstKE
LaunderSharmaKE
LienCubicKE
LienLeschziner
RNGkEpsilon
SSG
ShihQuadraticKE
SpalartAllmaras
kEpsilon
kOmega
kOmegaSST
kOmegaSSTLM
kOmegaSSTSAS
kkLOmega
qZeta
realizableKE
v2f
)
LES models
10
(
DeardorffDiffStress
Smagorinsky
SpalartAllmarasDDES
SpalartAllmarasDES
SpalartAllmarasIDDES
WALE
dynamicKEqn
dynamicLagrangian
kEqn
kOmegaSSTDES
)
Further work will be needed to support the -listTurbulenceModels option in
multiphase solvers.
103 lines
2.2 KiB
C
103 lines
2.2 KiB
C
bool listOptions = false ;
|
|
|
|
if
|
|
(
|
|
args.optionFound("listSwitches")
|
|
)
|
|
{
|
|
args.listSwitches(args.optionFound("includeUnsetSwitches"));
|
|
listOptions = true;
|
|
}
|
|
|
|
if
|
|
(
|
|
args.optionFound("listRegisteredSwitches")
|
|
)
|
|
{
|
|
args.listRegisteredSwitches(args.optionFound("includeUnsetSwitches"));
|
|
listOptions = true;
|
|
}
|
|
|
|
#ifdef fvPatchField_H
|
|
if (args.optionFound("listScalarBCs"))
|
|
{
|
|
Info<< "scalarBCs"
|
|
<< fvPatchField<scalar>::dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
listOptions = true;
|
|
}
|
|
|
|
if (args.optionFound("listVectorBCs"))
|
|
{
|
|
Info<< "vectorBCs"
|
|
<< fvPatchField<vector>::dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
listOptions = true;
|
|
}
|
|
#endif
|
|
|
|
#ifdef functionObject_H
|
|
if (args.optionFound("listFunctionObjects"))
|
|
{
|
|
Info<< "functionObjects"
|
|
<< functionObject::dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
listOptions = true;
|
|
}
|
|
#endif
|
|
|
|
#ifdef fvOption_H
|
|
if (args.optionFound("listFvOptions"))
|
|
{
|
|
Info<< "fvOptions"
|
|
<< fv::option::dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
listOptions = true;
|
|
}
|
|
#endif
|
|
|
|
#ifdef turbulentTransportModel_H
|
|
if (args.optionFound("listTurbulenceModels"))
|
|
{
|
|
Info<< "Turbulence models"
|
|
<< incompressible::turbulenceModel::
|
|
dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
|
|
Info<< "RAS models"
|
|
<< incompressible::RASModel::
|
|
dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
|
|
Info<< "LES models"
|
|
<< incompressible::LESModel::
|
|
dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
listOptions = true;
|
|
}
|
|
#elif defined(turbulentFluidThermoModel_H)
|
|
if (args.optionFound("listTurbulenceModels"))
|
|
{
|
|
Info<< "Turbulence models"
|
|
<< compressible::turbulenceModel::
|
|
dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
|
|
Info<< "RAS models"
|
|
<< compressible::RASModel::
|
|
dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
|
|
Info<< "LES models"
|
|
<< compressible::LESModel::
|
|
dictionaryConstructorTablePtr_->sortedToc()
|
|
<< endl;
|
|
listOptions = true;
|
|
}
|
|
#endif
|
|
|
|
if (listOptions)
|
|
{
|
|
exit(0);
|
|
}
|