ENH: consistent use of name scoping for turbulence fields

- with IOobject::scopedName instead of literal ':' for future
  transitioning of separators
This commit is contained in:
Mark Olesen 2023-05-03 19:55:40 +02:00
parent 8e32db2b5f
commit 2df90880d6
16 changed files with 62 additions and 64 deletions

View File

@ -7,7 +7,7 @@
-------------------------------------------------------------------------------
Copyright (C) 2011-2016 OpenFOAM Foundation
Copyright (C) 2022 Upstream CFD GmbH
Copyright (C) 2017-2022 OpenCFD Ltd.
Copyright (C) 2017-2023 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -131,15 +131,6 @@ class kOmegaSSTBase
:
public BasicEddyViscosityModel
{
// Private Member Functions
//- No copy construct
kOmegaSSTBase(const kOmegaSSTBase&) = delete;
//- No copy assignment
void operator=(const kOmegaSSTBase&) = delete;
protected:
// Protected Data
@ -239,7 +230,7 @@ protected:
{
return tmp<volScalarField::Internal>::New
(
this->type() + ":beta",
IOobject::scopedName(this->type(), "beta"),
blend(F1, beta1_, beta2_)
);
}
@ -251,7 +242,7 @@ protected:
{
return tmp<volScalarField::Internal>::New
(
this->type() + ":gamma",
IOobject::scopedName(this->type(), "gamma"),
blend(F1, gamma1_, gamma2_)
);
}
@ -313,6 +304,15 @@ public:
typedef typename BasicEddyViscosityModel::transportModel transportModel;
// Generated Methods
//- No copy construct
kOmegaSSTBase(const kOmegaSSTBase&) = delete;
//- No copy assignment
void operator=(const kOmegaSSTBase&) = delete;
// Constructors
//- Construct from components

View File

@ -264,7 +264,7 @@ void RNGkEpsilon<BasicTurbulenceModel>::correct()
tmp<volTensorField> tgradU = fvc::grad(U);
const volScalarField::Internal GbyNu
(
this->type() + ":GbyNu",
IOobject::scopedName(this->type(), "GbyNu"),
tgradU().v() && dev(twoSymm(tgradU().v()))
);
tgradU.clear();

View File

@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2017 OpenFOAM Foundation
Copyright (C) 2019-2020 OpenCFD Ltd.
Copyright (C) 2019-2023 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -245,7 +245,7 @@ void kEpsilon<BasicTurbulenceModel>::correct()
tmp<volTensorField> tgradU = fvc::grad(U);
const volScalarField::Internal GbyNu
(
this->type() + ":GbyNu",
IOobject::scopedName(this->type(), "GbyNu"),
tgradU().v() && dev(twoSymm(tgradU().v()))
);
const volScalarField::Internal G(this->GName(), nut()*GbyNu);

View File

@ -103,12 +103,11 @@ nu
const volScalarField& strainRate
) const
{
auto tnu =
volScalarField::New
(
IOobject::groupName(type() + ":nu", nu0.group()),
nu0.mesh(),
dimensionedScalar(dimViscosity, Zero)
auto tnu = volScalarField::New
(
IOobject::scopedName(type(), IOobject::groupName("nu", nu0.group())),
nu0.mesh(),
dimensionedScalar(dimViscosity, Zero)
);
tnu.ref().primitiveFieldRef() = strainRateFunction_->value(strainRate);

View File

@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2013-2015 OpenFOAM Foundation
Copyright (C) 2019-2020 OpenCFD Ltd.
Copyright (C) 2019-2023 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -38,8 +38,8 @@ SourceFiles
\*---------------------------------------------------------------------------*/
#ifndef turbulenceModel_H
#define turbulenceModel_H
#ifndef Foam_turbulenceModel_H
#define Foam_turbulenceModel_H
#include "IOdictionary.H"
#include "primitiveFieldsFwd.H"
@ -138,7 +138,7 @@ public:
//- Helper function to return the name of the turbulence G field
inline word GName() const
{
return word(type() + ":G");
return IOobject::scopedName(type(), "G");
}
//- Access function to velocity field

View File

@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2020 ENERCON GmbH
Copyright (C) 2020 OpenCFD Ltd
Copyright (C) 2020-2023 OpenCFD Ltd
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -71,7 +71,7 @@ void Foam::fv::atmAmbientTurbSource::atmAmbientTurbSourceOmega
const volScalarField::Internal& beta =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":beta")
IOobject::scopedName(turbPtr->type(), "beta")
);
// (RS:Eq. 4, rhs-term:5)

View File

@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2020 ENERCON GmbH
Copyright (C) 2020 OpenCFD Ltd.
Copyright (C) 2020-2023 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -52,7 +52,7 @@ void Foam::fv::atmBuoyancyTurbSource::atmBuoyancyTurbSourceEpsilon
const volScalarField::Internal& GbyNu =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":GbyNu")
IOobject::scopedName(turbPtr->type(), "GbyNu")
);
const volScalarField::Internal G(GbyNu*Cmu_*sqr(k())/epsilon());
@ -82,18 +82,18 @@ void Foam::fv::atmBuoyancyTurbSource::atmBuoyancyTurbSourceOmega
const volScalarField::Internal& GbyNu =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":GbyNu")
IOobject::scopedName(turbPtr->type(), "GbyNu")
);
const volScalarField::Internal G(GbyNu*Cmu_*k()/omega());
const volScalarField::Internal& gamma =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":gamma")
IOobject::scopedName(turbPtr->type(), "gamma")
);
const volScalarField::Internal& beta =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":beta")
IOobject::scopedName(turbPtr->type(), "beta")
);
// (ARAL:Eq. 5, rhs-term:3)

View File

@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2020 ENERCON GmbH
Copyright (C) 2020 OpenCFD Ltd
Copyright (C) 2020-2023 OpenCFD Ltd
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -52,7 +52,7 @@ void Foam::fv::atmLengthScaleTurbSource::atmLengthScaleTurbSourceEpsilon
const volScalarField::Internal& GbyNu =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":GbyNu")
IOobject::scopedName(turbPtr->type(), "GbyNu")
);
eqn += alpha()*rho()*calcC1Star(k, epsilon)*GbyNu*Cmu_*k;
@ -80,17 +80,17 @@ void Foam::fv::atmLengthScaleTurbSource::atmLengthScaleTurbSourceOmega
const volScalarField::Internal& GbyNu =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":GbyNu")
IOobject::scopedName(turbPtr->type(), "GbyNu")
);
const volScalarField::Internal& gamma =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":gamma")
IOobject::scopedName(turbPtr->type(), "gamma")
);
const volScalarField::Internal& beta =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":beta")
IOobject::scopedName(turbPtr->type(), "beta")
);
eqn += alpha()*rho()*calcGammaStar(k, omega, gamma, beta)*GbyNu;

View File

@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2020 ENERCON GmbH
Copyright (C) 2020 OpenCFD Ltd.
Copyright (C) 2020-2023 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -71,12 +71,12 @@ void Foam::fv::atmPlantCanopyTurbSource::atmPlantCanopyTurbSourceOmega
const volScalarField::Internal& gamma =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":gamma")
IOobject::scopedName(turbPtr->type(), "gamma")
);
const volScalarField::Internal& beta =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":beta")
IOobject::scopedName(turbPtr->type(), "beta")
);
eqn -= fvm::Sp(alpha()*rho()*(gamma - beta)*calcPlantCanopyTerm(U), omega);

View File

@ -5,7 +5,7 @@
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2020-2021 OpenCFD Ltd.
Copyright (C) 2020-2023 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -103,7 +103,7 @@ void Foam::fv::buoyancyTurbSource::buoyancyTurbSourceOmega
const auto& gamma =
mesh_.lookupObjectRef<volScalarField::Internal>
(
word(turbPtr->type() + ":gamma")
IOobject::scopedName(turbPtr->type(), "gamma")
);
// (Heuristic approximation to BMA:Eq. 6)

View File

@ -131,7 +131,7 @@ Foam::fv::interRegionHeatTransferModel::interRegionHeatTransferModel
(
IOobject
(
type() + ":htc",
IOobject::scopedName(type(), "htc"),
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
@ -169,19 +169,11 @@ void Foam::fv::interRegionHeatTransferModel::addSup
const auto& T = mesh_.lookupObject<volScalarField>(TName_);
auto tTmapped = tmp<volScalarField>::New
auto tTmapped = volScalarField::New
(
IOobject
(
type() + ":Tmapped",
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
IOobject::scopedName(type(), "Tmapped"),
T
);
auto& Tmapped = tTmapped.ref();
const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);

View File

@ -43,9 +43,14 @@ void Foam::SuppressionCollision<CloudType>::collide
lookupObject<kinematicCloud>(suppressionCloud_);
volScalarField vDotSweep(sc.vDotSweep());
dimensionedScalar Dt("dt", dimTime, dt);
volScalarField P(type() + ":p", 1.0 - exp(-vDotSweep*Dt));
auto tP = volScalarField::New
(
IOobject::scopedName(type(), "p"),
1.0 - exp(-vDotSweep*Dt)
);
const auto& P = tP();
for (typename CloudType::parcelType& p : this->owner())
{

View File

@ -1066,7 +1066,7 @@ void adjointkOmegaSST::updatePrimalRelatedFields()
(
IOobject
(
type() + ":G",
IOobject::scopedName(type(), "G"),
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,

View File

@ -50,8 +50,8 @@ SourceFiles
\*---------------------------------------------------------------------------*/
#ifndef adjointkOmegaSST_H
#define adjointkOmegaSST_H
#ifndef Foam_adjointkOmegaSST_H
#define Foam_adjointkOmegaSST_H
#include "adjointRASModel.H"
#include "wallDist.H"
@ -247,7 +247,7 @@ protected:
{
return tmp<volScalarField::Internal>::New
(
this->type() + ":beta",
IOobject::scopedName(this->type(), "beta"),
blend(F1, beta1_, beta2_)
);
}
@ -259,7 +259,7 @@ protected:
{
return tmp<volScalarField>::New
(
this->type() + ":beta",
IOobject::scopedName(this->type(), "beta"),
blend(F1, beta1_, beta2_)
);
}
@ -271,7 +271,7 @@ protected:
{
return tmp<volScalarField::Internal>::New
(
this->type() + ":gamma",
IOobject::scopedName(this->type(), "gamma"),
blend(F1, gamma1_, gamma2_)
);
}
@ -283,7 +283,7 @@ protected:
{
return tmp<volScalarField>::New
(
this->type() + ":gamma",
IOobject::scopedName(this->type(), "gamma"),
blend(F1, gamma1_, gamma2_)
);
}

View File

@ -113,9 +113,11 @@ tmp<volScalarField::Internal> kOmegaSST::computeG()
tmp<volTensorField> tgradU = fvc::grad(U);
volScalarField::Internal GbyNu0
(
this->type() + ":GbyNu",
IOobject::scopedName(this->type(), "GbyNu"),
(tgradU() && dev(twoSymm(tgradU())))
);
// NB: leave tmp registered (for correctBoundaryConditions)
auto tG =
tmp<volScalarField::Internal>::New
(

View File

@ -72,7 +72,7 @@ Foam::diameterModels::driftModels::phaseChange::phaseChange
(
IOobject
(
IOobject::groupName(type() + ":W", pair.name()),
IOobject::scopedName(type(), IOobject::groupName("W", pair.name())),
popBal_.mesh().time().timeName(),
popBal_.mesh()
),