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-------------------------------------------------------------------------------
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
29#include "phaseSystem.H"
31#include "dictionary.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace functionObjects
38{
41 (
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
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 {
87 const volScalarField& he = phase.thermo().he();
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 =
144 (
146 (
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// ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
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.
Definition: Field.H:82
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
word resultName_
Name of result field.
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 reactingEuler solvers.
virtual bool calc()
Calculate the heat transfer coefficient field.
virtual bool read(const dictionary &dict)
Read the heatTransferCoeff data.
tmp< FieldField< Field, scalar > > q() const
Calculate heat flux.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
const Type & lookupObject(const word &name, const bool recursive=false) const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:76
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:37
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
dynamicFvMesh & mesh
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Namespace for OpenFOAM.
const dimensionSet dimPower
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & alpha
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
Definition: createFields.H:66
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333