mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: prefer 'model' keyword instead of 'RASModel',... (#149)
- since the context (laminar/RAS/LES) is already given by the sub-dictionary, it is redundant to use as prefix as well. - silently support the longer names as compat methods
This commit is contained in:
@ -6,7 +6,7 @@
|
|||||||
\\/ M anipulation |
|
\\/ M anipulation |
|
||||||
-------------------------------------------------------------------------------
|
-------------------------------------------------------------------------------
|
||||||
Copyright (C) 2015 OpenFOAM Foundation
|
Copyright (C) 2015 OpenFOAM Foundation
|
||||||
Copyright (C) 2015-2016 OpenCFD Ltd.
|
Copyright (C) 2015-2020 OpenCFD Ltd.
|
||||||
-------------------------------------------------------------------------------
|
-------------------------------------------------------------------------------
|
||||||
License
|
License
|
||||||
This file is part of OpenFOAM.
|
This file is part of OpenFOAM.
|
||||||
@ -93,17 +93,18 @@ Foam::dictionary Foam::solverTemplate::readFluidFieldTemplates
|
|||||||
|
|
||||||
const dictionary fieldModels(solverDict.subDict("fluidModels"));
|
const dictionary fieldModels(solverDict.subDict("fluidModels"));
|
||||||
|
|
||||||
word turbulenceModel("laminar"); // default to laminar
|
word turbModel("laminar"); // default to laminar
|
||||||
word turbulenceType("none");
|
word turbType("none");
|
||||||
|
|
||||||
if (fieldModels.readIfPresent("turbulenceModel", turbulenceType))
|
if (fieldModels.readIfPresent("turbulenceModel", turbType))
|
||||||
{
|
{
|
||||||
if (turbulenceType == "turbulenceModel")
|
if (turbType == "turbulenceModel")
|
||||||
{
|
{
|
||||||
IOdictionary turbPropDict
|
IOdictionary turbPropDict
|
||||||
(
|
(
|
||||||
IOobject
|
IOobject
|
||||||
(
|
(
|
||||||
|
// turbulenceModel::propertiesName
|
||||||
"turbulenceProperties",
|
"turbulenceProperties",
|
||||||
runTime.constant(),
|
runTime.constant(),
|
||||||
regionName,
|
regionName,
|
||||||
@ -118,17 +119,19 @@ Foam::dictionary Foam::solverTemplate::readFluidFieldTemplates
|
|||||||
|
|
||||||
if (modelType == "laminar")
|
if (modelType == "laminar")
|
||||||
{
|
{
|
||||||
// Leave turbulenceModel as laminar
|
// Leave as laminar
|
||||||
}
|
}
|
||||||
else if (modelType == "RAS")
|
else if (modelType == "RAS")
|
||||||
{
|
{
|
||||||
|
// "RASModel" for v2006 and earlier
|
||||||
turbPropDict.subDict(modelType)
|
turbPropDict.subDict(modelType)
|
||||||
.readEntry("RASModel", turbulenceModel);
|
.readCompat("model", {{"RASModel", -2006}}, turbModel);
|
||||||
}
|
}
|
||||||
else if (modelType == "LES")
|
else if (modelType == "LES")
|
||||||
{
|
{
|
||||||
|
// "LESModel" for v2006 and earlier
|
||||||
turbPropDict.subDict(modelType)
|
turbPropDict.subDict(modelType)
|
||||||
.readEntry("LESModel", turbulenceModel);
|
.readCompat("model", {{"LESModel", -2006}}, turbModel);
|
||||||
}
|
}
|
||||||
else
|
else
|
||||||
{
|
{
|
||||||
@ -141,20 +144,19 @@ Foam::dictionary Foam::solverTemplate::readFluidFieldTemplates
|
|||||||
else
|
else
|
||||||
{
|
{
|
||||||
FatalErrorInFunction
|
FatalErrorInFunction
|
||||||
<< "Unhandled turbulence model option " << turbulenceType
|
<< "Unhandled turbulence model option " << turbType
|
||||||
<< ". Valid options are turbulenceModel"
|
<< ". Valid options are turbulenceModel"
|
||||||
<< exit(FatalError);
|
<< exit(FatalError);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
Info<< " Selecting " << turbulenceType << ": " << turbulenceModel
|
Info<< " Selecting " << turbType << ": " << turbModel << endl;
|
||||||
<< endl;
|
|
||||||
|
|
||||||
IOdictionary turbModelDict
|
IOdictionary turbModelDict
|
||||||
(
|
(
|
||||||
IOobject
|
IOobject
|
||||||
(
|
(
|
||||||
fileName(turbModelDir/turbulenceModel),
|
fileName(turbModelDir/turbModel),
|
||||||
runTime,
|
runTime,
|
||||||
IOobject::MUST_READ
|
IOobject::MUST_READ
|
||||||
)
|
)
|
||||||
|
|||||||
@ -151,7 +151,11 @@ Foam::LESModel<BasicTurbulenceModel>::New
|
|||||||
|
|
||||||
const dictionary& dict = modelDict.subDict("LES");
|
const dictionary& dict = modelDict.subDict("LES");
|
||||||
|
|
||||||
const word modelType(dict.get<word>("LESModel"));
|
const word modelType
|
||||||
|
(
|
||||||
|
// "LESModel" for v2006 and earlier
|
||||||
|
dict.getCompat<word>("model", {{"LESModel", -2006}})
|
||||||
|
);
|
||||||
|
|
||||||
Info<< "Selecting LES turbulence model " << modelType << endl;
|
Info<< "Selecting LES turbulence model " << modelType << endl;
|
||||||
|
|
||||||
@ -162,7 +166,7 @@ Foam::LESModel<BasicTurbulenceModel>::New
|
|||||||
FatalIOErrorInLookup
|
FatalIOErrorInLookup
|
||||||
(
|
(
|
||||||
dict,
|
dict,
|
||||||
"LESModel",
|
"LES model",
|
||||||
modelType,
|
modelType,
|
||||||
*dictionaryConstructorTablePtr_
|
*dictionaryConstructorTablePtr_
|
||||||
) << exit(FatalIOError);
|
) << exit(FatalIOError);
|
||||||
|
|||||||
@ -141,7 +141,11 @@ Foam::RASModel<BasicTurbulenceModel>::New
|
|||||||
|
|
||||||
const dictionary& dict = modelDict.subDict("RAS");
|
const dictionary& dict = modelDict.subDict("RAS");
|
||||||
|
|
||||||
const word modelType(dict.get<word>("RASModel"));
|
const word modelType
|
||||||
|
(
|
||||||
|
// "RASModel" for v2006 and earlier
|
||||||
|
dict.getCompat<word>("model", {{"RASModel", -2006}})
|
||||||
|
);
|
||||||
|
|
||||||
Info<< "Selecting RAS turbulence model " << modelType << endl;
|
Info<< "Selecting RAS turbulence model " << modelType << endl;
|
||||||
|
|
||||||
@ -152,7 +156,7 @@ Foam::RASModel<BasicTurbulenceModel>::New
|
|||||||
FatalIOErrorInLookup
|
FatalIOErrorInLookup
|
||||||
(
|
(
|
||||||
dict,
|
dict,
|
||||||
"RASModel",
|
"RAS model",
|
||||||
modelType,
|
modelType,
|
||||||
*dictionaryConstructorTablePtr_
|
*dictionaryConstructorTablePtr_
|
||||||
) << exit(FatalIOError);
|
) << exit(FatalIOError);
|
||||||
|
|||||||
@ -112,7 +112,11 @@ Foam::laminarModel<BasicTurbulenceModel>::New
|
|||||||
{
|
{
|
||||||
const dictionary& dict = *dictptr;
|
const dictionary& dict = *dictptr;
|
||||||
|
|
||||||
const word modelType(dict.get<word>("laminarModel"));
|
const word modelType
|
||||||
|
(
|
||||||
|
// laminarModel -> model (after v2006)
|
||||||
|
dict.getCompat<word>("model", {{"laminarModel", -2006}})
|
||||||
|
);
|
||||||
|
|
||||||
Info<< "Selecting laminar stress model " << modelType << endl;
|
Info<< "Selecting laminar stress model " << modelType << endl;
|
||||||
|
|
||||||
@ -123,7 +127,7 @@ Foam::laminarModel<BasicTurbulenceModel>::New
|
|||||||
FatalIOErrorInLookup
|
FatalIOErrorInLookup
|
||||||
(
|
(
|
||||||
dict,
|
dict,
|
||||||
"laminarModel",
|
"laminar model",
|
||||||
modelType,
|
modelType,
|
||||||
*dictionaryConstructorTablePtr_
|
*dictionaryConstructorTablePtr_
|
||||||
) << exit(FatalIOError);
|
) << exit(FatalIOError);
|
||||||
|
|||||||
@ -272,9 +272,19 @@ autoPtr<RASModelVariables> RASModelVariables::New
|
|||||||
)
|
)
|
||||||
);
|
);
|
||||||
|
|
||||||
const dictionary dict(modelDict.subOrEmptyDict("RAS"));
|
word modelType("laminar"); // default to laminar
|
||||||
|
|
||||||
const word modelType(dict.getOrDefault<word>("RASModel", "laminar"));
|
const dictionary* dictptr = modelDict.findDict("RAS");
|
||||||
|
|
||||||
|
if (dictptr)
|
||||||
|
{
|
||||||
|
// "RASModel" for v2006 and earlier
|
||||||
|
dictptr->readCompat("model", {{"RASModel", -2006}}, modelType);
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
dictptr = &dictionary::null;
|
||||||
|
}
|
||||||
|
|
||||||
Info<< "Creating references for RASModel variables : " << modelType << endl;
|
Info<< "Creating references for RASModel variables : " << modelType << endl;
|
||||||
|
|
||||||
@ -284,7 +294,7 @@ autoPtr<RASModelVariables> RASModelVariables::New
|
|||||||
{
|
{
|
||||||
FatalIOErrorInLookup
|
FatalIOErrorInLookup
|
||||||
(
|
(
|
||||||
dict,
|
*dictptr,
|
||||||
"RASModelVariables",
|
"RASModelVariables",
|
||||||
modelType,
|
modelType,
|
||||||
*dictionaryConstructorTablePtr_
|
*dictionaryConstructorTablePtr_
|
||||||
|
|||||||
Reference in New Issue
Block a user