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-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
29 #include "heatTransferCoeffModel.H"
30 #include "multiphaseInterSystem.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(multiphaseInterHtcModel, 0);
41  (
42  functionObject,
43  multiphaseInterHtcModel,
44  dictionary
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 
72  auto tq = tmp<FieldField<Field, scalar>>::New(Tbf.size());
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  }
111 
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  setResultName(typeName, "htc:" + htcModelPtr_->type());
131 
132  auto* htcPtr =
133  new volScalarField
134  (
135  IOobject
136  (
137  resultName_,
138  mesh_.time().timeName(),
139  mesh_,
143  ),
144  mesh_,
146  );
147 
148  mesh_.objectRegistry::store(htcPtr);
149 }
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
155 (
156  const dictionary& dict
157 )
158 {
159  if (!fieldExpression::read(dict) || !htcModelPtr_->read(dict))
160  {
161  return false;
162  }
163 
164  return true;
165 }
166 
167 
168 // ************************************************************************* //
twoPhaseSystem & fluid
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
word resultName_
Name of result field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual bool calc()
Calculate the heat transfer coefficient field.
engineTime & runTime
virtual bool read(const dictionary &dict)
Read the fieldExpression data.
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.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
A base class for heat transfer coefficient models.
tmp< FieldField< Field, scalar > > q() const
Calculate heat flux.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
virtual bool read(const dictionary &dict)
Read the heatTransferCoeff data.
multiphaseInterHtcModel(const multiphaseInterHtcModel &)=delete
No copy construct.
const dimensionSet dimPower
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
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field T\"<< endl;volScalarField T(IOobject("T", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Calculating field g.h\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), p_rgh);Info<< "Creating multiphaseSystem\"<< endl;autoPtr< multiphaseInter::multiphaseSystem > fluidPtr
Definition: createFields.H:66
Nothing to be read.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
const fvMesh & mesh_
Reference to the fvMesh.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127