Loading...
Searching...
No Matches
multiphaseInterHtcModel.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2022-2024 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace functionObjects
38{
41 (
45 );
46}
47}
48
49
50// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51
53{
54 auto& htc =
55 htcModelPtr_->mesh().lookupObjectRef<volScalarField>(resultName_);
56
57 htcModelPtr_->calc(htc, q());
58
59 return true;
60}
61
62
65{
66 const fvMesh& mesh = htcModelPtr_->mesh();
67
68 const auto& T = mesh.lookupObject<volScalarField>(htcModelPtr_->TName());
69
70 const volScalarField::Boundary& Tbf = T.boundaryField();
71
73 auto& q = tq.ref();
74
75 forAll(q, patchi)
76 {
77 q.set(patchi, new Field<scalar>(Tbf[patchi].size(), Zero));
78 }
79
80 const auto* fluidPtr =
81 mesh.cfindObject<multiphaseInterSystem>("phaseProperties");
82
83 if (!fluidPtr)
84 {
86 << "Unable to find a valid phaseSystem to evaluate q" << nl
87 << exit(FatalError);
88 }
89
90 const multiphaseInterSystem& fluid = *fluidPtr;
91
92 for (const label patchi : htcModelPtr_->patchIDs())
93 {
94 q[patchi] += fluid.kappaEff(patchi)()*Tbf[patchi].snGrad();
95 }
96
97 // Add radiative heat flux contribution if present
98
99 const auto* qrPtr =
100 mesh.cfindObject<volScalarField>(htcModelPtr_->qrName());
101
102 if (qrPtr)
103 {
104 const volScalarField::Boundary& qrbf = qrPtr->boundaryField();
105
106 for (const label patchi : htcModelPtr_->patchIDs())
107 {
108 q[patchi] += qrbf[patchi];
109 }
110 }
112 return tq;
113}
114
115
116// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117
119(
120 const word& name,
121 const Time& runTime,
122 const dictionary& dict
123)
124:
126 htcModelPtr_(heatTransferCoeffModel::New(dict, mesh_, fieldName_))
127{
128 read(dict);
129
130 // scopedName?
131 setResultName(typeName, "htc:" + htcModelPtr_->type());
132
133 auto* htcPtr =
135 (
137 (
139 mesh_.time().timeName(),
140 mesh_.thisDb(),
144 ),
145 mesh_,
147 );
149 regIOobject::store(htcPtr);
150}
151
152
153// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154
156(
157 const dictionary& dict
158)
159{
160 if (!fieldExpression::read(dict) || !htcModelPtr_->read(dict))
161 {
162 return false;
163 }
164
165 return true;
166}
167
168
169// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
twoPhaseSystem & fluid
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
word resultName_
Name of result field.
word fieldName_
Name of field to process.
fieldExpression(const word &name, const Time &runTime, const dictionary &dict, const word &fieldName=word::null, const word &resultName=word::null)
Construct from name, Time and dictionary.
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
const fvMesh & mesh_
Reference to the fvMesh.
A heat transfer coefficient for multiphase inter solvers (e.g. icoReactingMultiphaseFoam).
virtual bool calc()
Calculate the heat transfer coefficient field.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
multiphaseInterHtcModel(const multiphaseInterHtcModel &)=delete
No copy construct.
tmp< FieldField< Field, scalar > > q() const
Calculate heat flux.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
A base class for heat transfer coefficient models.
bool store()
Register object with its registry and transfer ownership to the registry.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
const dimensionSet dimPower
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
Info<< "Reading field p_rgh\n"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field T\n"<< endl;volScalarField T(IOobject("T", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Calculating field g.h\n"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), p_rgh);Info<< "Creating multiphaseSystem\n"<< endl;autoPtr< multiphaseInter::multiphaseSystem > fluidPtr
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299