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 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 "phaseSystem.H"
31 #include "dictionary.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
41  (
42  functionObject,
44  dictionary
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51 
54 {
55  const fvMesh& mesh = htcModelPtr_->mesh();
56 
57  const volScalarField& T =
58  mesh.lookupObject<volScalarField>(htcModelPtr_->TName());
59 
60  const volScalarField::Boundary& Tbf = T.boundaryField();
61 
62  auto tq = tmp<FieldField<Field, scalar>>::New(Tbf.size());
63  auto& q = tq.ref();
64 
65  forAll(q, patchi)
66  {
67  q.set(patchi, new Field<scalar>(Tbf[patchi].size(), Zero));
68  }
69 
70  const auto* fluidPtr =
71  mesh.cfindObject<phaseSystem>("phaseProperties");
72 
73  if (!fluidPtr)
74  {
76  << "Unable to find a valid phaseSystem to evaluate q" << nl
77  << exit(FatalError);
78  }
79 
80  const phaseSystem& fluid = *fluidPtr;
81 
82  for (const label patchi : htcModelPtr_->patchSet())
83  {
84  for (const phaseModel& phase : fluid.phases())
85  {
86  const fvPatchScalarField& alpha = phase.boundaryField()[patchi];
87  const volScalarField& he = phase.thermo().he();
88  const volScalarField::Boundary& hebf = he.boundaryField();
89 
90  q[patchi] +=
91  alpha*phase.alphaEff(patchi)()*hebf[patchi].snGrad();
92  }
93  }
94 
95  // Add radiative heat flux contribution if present
96 
97  const volScalarField* qrPtr =
98  mesh.cfindObject<volScalarField>(htcModelPtr_->qrName());
99 
100  if (qrPtr)
101  {
102  const volScalarField::Boundary& qrbf = qrPtr->boundaryField();
103 
104  for (const label patchi : htcModelPtr_->patchSet())
105  {
106  q[patchi] += qrbf[patchi];
107  }
108  }
109 
110  return tq;
111 }
112 
113 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
114 
116 {
117  auto& htc =
118  htcModelPtr_->mesh().lookupObjectRef<volScalarField>(resultName_);
119 
120  htcModelPtr_->calc(htc, q());
121 
122  return true;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
129 (
130  const word& name,
131  const Time& runTime,
132  const dictionary& dict
133 )
134 :
136  htcModelPtr_(nullptr)
137 {
138  read(dict);
139 
140  setResultName(typeName, "htc:" + htcModelPtr_->type());
141 
142  volScalarField* htcPtr =
143  new volScalarField
144  (
145  IOobject
146  (
147  resultName_,
148  mesh_.time().timeName(),
149  mesh_,
152  ),
153  mesh_,
155  );
156 
157  mesh_.objectRegistry::store(htcPtr);
158 }
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
164 {
166  {
167  htcModelPtr_ = heatTransferCoeffModel::New(dict, mesh_, fieldName_);
168 
169  htcModelPtr_->read(dict);
170 
171  return true;
172  }
173 
174  return false;
175 }
176 
177 
178 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::functionObjects::reactingEulerHtcModel::q
tmp< FieldField< Field, scalar > > q() const
Calculate heat flux.
Definition: reactingEulerHtcModel.C:53
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
reactingEulerHtcModel
A heat transfer coefficient for reactingEuler solvers.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::functionObject::New
static autoPtr< functionObject > New(const word &name, const Time &runTime, const dictionary &dict)
Select from dictionary, based on its "type" entry.
Definition: functionObject.C:79
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::functionObjects::fieldExpression::read
virtual bool read(const dictionary &dict)
Read the fieldExpression data.
Definition: fieldExpression.C:93
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
fluidPtr
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< multiphaseSystem > fluidPtr
Definition: createFields.H:66
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::Field< scalar >
reactingEulerHtcModel.H
Foam::heatTransferCoeffModel::New
static autoPtr< heatTransferCoeffModel > New(const dictionary &dict, const fvMesh &mesh, const word &TName)
Return a reference to the selected heat transfer coefficient model.
Definition: heatTransferCoeffModelNew.C:35
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::dimPower
const dimensionSet dimPower
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::functionObjects::reactingEulerHtcModel::calc
virtual bool calc()
Calculate the heat transfer coefficient field.
Definition: reactingEulerHtcModel.C:115
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::functionObjects::fieldExpression
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
Definition: fieldExpression.H:120
he
volScalarField & he
Definition: YEEqn.H:52
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
dictionary.H
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::functionObjects::reactingEulerHtcModel::reactingEulerHtcModel
reactingEulerHtcModel(const reactingEulerHtcModel &)=delete
No copy construct.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:66
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::functionObjects::reactingEulerHtcModel::read
virtual bool read(const dictionary &dict)
Read the heatTransferCoeff data.
Definition: reactingEulerHtcModel.C:163
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62