From dbdae526933e533dd71694aeccfa4d8de1956c51 Mon Sep 17 00:00:00 2001 From: Kutalmis Bercin <kutalmis.bercin@esi-group.com> Date: Fri, 25 Jun 2021 19:08:55 +0100 Subject: [PATCH] INT: turbulenceFields: add phase support ENH: turbulenceFields: add const specifier to compressible func --- .../field/turbulenceFields/turbulenceFields.C | 38 ++++++++++++++----- .../field/turbulenceFields/turbulenceFields.H | 9 ++++- .../turbulenceFieldsTemplates.C | 7 +++- 3 files changed, 40 insertions(+), 14 deletions(-) diff --git a/src/functionObjects/field/turbulenceFields/turbulenceFields.C b/src/functionObjects/field/turbulenceFields/turbulenceFields.C index 894f8adf8f8..b01ae3695ac 100644 --- a/src/functionObjects/field/turbulenceFields/turbulenceFields.C +++ b/src/functionObjects/field/turbulenceFields/turbulenceFields.C @@ -5,7 +5,7 @@ \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- - Copyright (C) 2013-2016 OpenFOAM Foundation + Copyright (C) 2013-2020 OpenFOAM Foundation Copyright (C) 2015-2021 OpenCFD Ltd. ------------------------------------------------------------------------------- License @@ -94,7 +94,10 @@ void Foam::functionObjects::turbulenceFields::initialise() { for (const word& f : fieldSet_) { - const word scopedName(prefix_ + f); + const word scopedName + ( + IOobject::groupName(prefix_ + f, phaseName_) + ); if (obr_.found(scopedName)) { @@ -111,13 +114,16 @@ void Foam::functionObjects::turbulenceFields::initialise() } -bool Foam::functionObjects::turbulenceFields::compressible() +bool Foam::functionObjects::turbulenceFields::compressible +( + const word& modelName +) const { - if (obr_.foundObject<compressible::turbulenceModel>(modelName_)) + if (obr_.foundObject<compressible::turbulenceModel>(modelName)) { return true; } - else if (obr_.foundObject<incompressible::turbulenceModel>(modelName_)) + else if (obr_.foundObject<incompressible::turbulenceModel>(modelName)) { return false; } @@ -142,6 +148,7 @@ Foam::functionObjects::turbulenceFields::turbulenceFields fvMeshFunctionObject(name, runTime, dict), initialised_(false), prefix_(dict.getOrDefault<word>("prefix", "turbulenceProperties:")), + phaseName_(dict.getOrDefault<word>("phase", word::null)), fieldSet_() { read(dict); @@ -155,6 +162,7 @@ bool Foam::functionObjects::turbulenceFields::read(const dictionary& dict) if (fvMeshFunctionObject::read(dict)) { dict.readIfPresent("prefix", prefix_); + dict.readIfPresent("phase", phaseName_); if (dict.found("field")) { @@ -171,7 +179,9 @@ bool Foam::functionObjects::turbulenceFields::read(const dictionary& dict) Info<< "storing fields:" << nl; for (const word& f : fieldSet_) { - Info<< " " << prefix_ << f << nl; + Info<< " " + << IOobject::groupName(prefix_ + f, phaseName_) + << nl; } Info<< endl; } @@ -196,12 +206,17 @@ bool Foam::functionObjects::turbulenceFields::execute() initialise(); } - const bool comp = compressible(); + const word modelName + ( + IOobject::groupName(modelName_, phaseName_) + ); + + const bool comp = compressible(modelName); if (comp) { const auto& model = - obr_.lookupObject<compressible::turbulenceModel>(modelName_); + obr_.lookupObject<compressible::turbulenceModel>(modelName); for (const word& f : fieldSet_) { @@ -278,7 +293,7 @@ bool Foam::functionObjects::turbulenceFields::execute() else { const auto& model = - obr_.lookupObject<incompressible::turbulenceModel>(modelName_); + obr_.lookupObject<incompressible::turbulenceModel>(modelName); for (const word& f : fieldSet_) { @@ -351,7 +366,10 @@ bool Foam::functionObjects::turbulenceFields::write() { for (const word& f : fieldSet_) { - const word scopedName(prefix_ + f); + const word scopedName + ( + IOobject::groupName(prefix_ + f, phaseName_) + ); writeObject(scopedName); } diff --git a/src/functionObjects/field/turbulenceFields/turbulenceFields.H b/src/functionObjects/field/turbulenceFields/turbulenceFields.H index 05ed772e476..f853d52e1e5 100644 --- a/src/functionObjects/field/turbulenceFields/turbulenceFields.H +++ b/src/functionObjects/field/turbulenceFields/turbulenceFields.H @@ -5,7 +5,7 @@ \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- - Copyright (C) 2013-2016 OpenFOAM Foundation + Copyright (C) 2013-2020 OpenFOAM Foundation Copyright (C) 2015-2021 OpenCFD Ltd. ------------------------------------------------------------------------------- License @@ -80,6 +80,7 @@ Usage // Optional entries (runtime modifiable) prefix <word>; + phase <word>; // Inherited entries ... @@ -94,6 +95,7 @@ Usage fields | Names of fields to store (see below) | wordList | yes | - field | Name of a field to store (see below) | word | yes | - prefix | Name of output-field prefix | word | no | turbulenceProperties: + phase | Name of operand phase | word | no | - \endtable where \c fields can include: @@ -214,6 +216,9 @@ protected: //- Name of output-field prefix word prefix_; + //- Name of operand phase + word phaseName_; + //- Fields to load wordHashSet fieldSet_; @@ -224,7 +229,7 @@ protected: void initialise(); //- Return true if compressible turbulence model is identified - bool compressible(); + bool compressible(const word& modelName) const; //- Process the turbulence field template<class Type> diff --git a/src/functionObjects/field/turbulenceFields/turbulenceFieldsTemplates.C b/src/functionObjects/field/turbulenceFields/turbulenceFieldsTemplates.C index a5452640d57..5f17e494343 100644 --- a/src/functionObjects/field/turbulenceFields/turbulenceFieldsTemplates.C +++ b/src/functionObjects/field/turbulenceFields/turbulenceFieldsTemplates.C @@ -5,7 +5,7 @@ \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- - Copyright (C) 2012-2016 OpenFOAM Foundation + Copyright (C) 2012-2020 OpenFOAM Foundation Copyright (C) 2018-2021 OpenCFD Ltd. ------------------------------------------------------------------------------- License @@ -40,7 +40,10 @@ void Foam::functionObjects::turbulenceFields::processField { typedef GeometricField<Type, fvPatchField, volMesh> FieldType; - const word scopedName(prefix_ + fieldName); + const word scopedName + ( + IOobject::groupName(prefix_ + fieldName, phaseName_) + ); FieldType* fldPtr = obr_.getObjectPtr<FieldType>(scopedName); -- GitLab