reactingEulerHtcModel.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) 2020-2022 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 
28 #include "reactingEulerHtcModel.H"
29 #include "heatTransferCoeffModel.H"
30 #include "phaseSystem.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(reactingEulerHtcModel, 0);
41  (
42  functionObject,
43  reactingEulerHtcModel,
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 volScalarField& T =
69  mesh.lookupObject<volScalarField>(htcModelPtr_->TName());
70 
71  const volScalarField::Boundary& Tbf = T.boundaryField();
72 
73  auto tq = tmp<FieldField<Field, scalar>>::New(Tbf.size());
74  auto& q = tq.ref();
75 
76  forAll(q, patchi)
77  {
78  q.set(patchi, new Field<scalar>(Tbf[patchi].size(), Zero));
79  }
80 
81  const auto* fluidPtr =
82  mesh.cfindObject<phaseSystem>("phaseProperties");
83 
84  if (!fluidPtr)
85  {
87  << "Unable to find a valid phaseSystem to evaluate q" << nl
88  << exit(FatalError);
89  }
90 
91  const phaseSystem& fluid = *fluidPtr;
92 
93  for (const label patchi : htcModelPtr_->patchIDs())
94  {
95  for (const phaseModel& phase : fluid.phases())
96  {
97  const fvPatchScalarField& alpha = phase.boundaryField()[patchi];
98  const volScalarField& he = phase.thermo().he();
99  const volScalarField::Boundary& hebf = he.boundaryField();
100 
101  q[patchi] +=
102  alpha*phase.alphaEff(patchi)()*hebf[patchi].snGrad();
103  }
104  }
105 
106  // Add radiative heat flux contribution if present
107 
108  const volScalarField* qrPtr =
109  mesh.cfindObject<volScalarField>(htcModelPtr_->qrName());
110 
111  if (qrPtr)
112  {
113  const volScalarField::Boundary& qrbf = qrPtr->boundaryField();
114 
115  for (const label patchi : htcModelPtr_->patchIDs())
116  {
117  q[patchi] += qrbf[patchi];
118  }
119  }
120 
121  return tq;
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
126 
128 (
129  const word& name,
130  const Time& runTime,
131  const dictionary& dict
132 )
133 :
135  htcModelPtr_(heatTransferCoeffModel::New(dict, mesh_, fieldName_))
136 {
137  read(dict);
138 
139  setResultName(typeName, "htc:" + htcModelPtr_->type());
140 
141  volScalarField* htcPtr =
142  new volScalarField
143  (
144  IOobject
145  (
146  resultName_,
147  mesh_.time().timeName(),
148  mesh_,
152  ),
153  mesh_,
155  );
157  mesh_.objectRegistry::store(htcPtr);
158 }
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
164 {
165  if (!fieldExpression::read(dict) || htcModelPtr_->read(dict))
166  {
167  return false;
168  }
169 
170  return true;
171 }
172 
173 
174 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
twoPhaseSystem & fluid
reactingEulerHtcModel(const reactingEulerHtcModel &)=delete
No copy construct.
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
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
virtual bool calc()
Calculate the heat transfer coefficient field.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
fvPatchField< scalar > fvPatchScalarField
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
virtual bool read(const dictionary &dict)
Read the heatTransferCoeff data.
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
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.
tmp< FieldField< Field, scalar > > q() const
Calculate heat flux.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
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