uniformFixedMultiphaseHeatFluxFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "fvPatchFieldMapper.H"
28 #include "phaseSystem.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const fvPatch& p,
38  const dictionary& dict
39 )
40 :
41  mixedFvPatchScalarField(p, iF, dict, false),
42  q_(Function1<scalar>::New("q", dict)),
43  relax_(dict.lookupOrDefault<scalar>("relax", 1))
44 {
45  valueFraction() = 1;
46  refValue() = patchInternalField();
47  refGrad() = Zero;
48 
49  operator==(patchInternalField());
50 }
51 
52 
55 (
57  const fvPatch& p,
59  const fvPatchFieldMapper& mapper
60 )
61 :
62  mixedFvPatchScalarField(psf, p, iF, mapper),
63  q_(psf.q_, false),
64  relax_(psf.relax_)
65 {}
66 
67 
70 (
73 )
74 :
75  mixedFvPatchScalarField(psf, iF),
76  q_(psf.q_, false),
77  relax_(psf.relax_)
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
84 {
85  if (updated())
86  {
87  return;
88  }
89 
90  const label patchi = patch().index();
91 
92  const scalar q = q_->value(db().time().userTimeValue());
93 
94  const phaseSystem& fluid =
95  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
96 
97  const phaseModel& thisPhase = fluid.phases()[internalField().group()];
98 
99  // Sums of alpha*kappaEff
100  scalarField sumAlphaKappaEff(patch().size(), rootVSmall);
101  scalarField sumNotThisAlphaKappaEff(patch().size(), rootVSmall);
102  scalarField sumNotThisAlphaKappaEffT(patch().size(), rootVSmall);
103 
104  // Contributions from phases other than this one
105  forAll(fluid.phases(), phasei)
106  {
107  const phaseModel& phase = fluid.phases()[phasei];
108 
109  if (&phase != &thisPhase)
110  {
111  const scalarField& alpha = phase.boundaryField()[patchi];
112  const scalarField kappaEff(phase.kappaEff(patchi));
113  const scalarField alphaKappaEff(alpha*kappaEff);
114  const fvPatchScalarField& T =
115  phase.thermo().T().boundaryField()[patchi];
116 
117  sumAlphaKappaEff += alphaKappaEff;
118  sumNotThisAlphaKappaEff += alphaKappaEff;
119  sumNotThisAlphaKappaEffT += alphaKappaEff*T.patchInternalField();
120  }
121  }
122 
123  // Contribution from this phase, and stabilisation
124  const scalarField& alpha = thisPhase.boundaryField()[patchi];
125  const scalarField kappaEff(thisPhase.kappaEff(patchi));
126  const scalarField alphaKappaEff(alpha*kappaEff);
127  const fvPatchScalarField& T =
128  thisPhase.thermo().T().boundaryField()[patchi];
129 
130  sumAlphaKappaEff += alphaKappaEff;
131  sumNotThisAlphaKappaEff =
132  max
133  (
134  sumNotThisAlphaKappaEff,
135  rootSmall*kappaEff
136  );
137  sumNotThisAlphaKappaEffT =
138  max
139  (
140  sumNotThisAlphaKappaEffT,
141  rootSmall*kappaEff*T.patchInternalField()
142  );
143 
144  // Mixed parameters
145  valueFraction() = sumNotThisAlphaKappaEff/sumAlphaKappaEff;
146  refValue() = sumNotThisAlphaKappaEffT/sumNotThisAlphaKappaEff;
147  refGrad() = q/max(alpha, rootSmall)/kappaEff;
148 
149  // Modify mixed parameters for under-relaxation
150  if (relax_ != 1)
151  {
152  const scalarField f(valueFraction());
153  valueFraction() = 1 - relax_*(1 - f);
154  refValue() = (f*relax_*refValue() + (1 - relax_)*T)/valueFraction();
155  //refGrad() = refGrad(); // No change
156  }
157 
158  mixedFvPatchScalarField::updateCoeffs();
159 }
160 
161 
163 (
164  Ostream& os
165 ) const
166 {
168  writeEntry(os, q_());
169  writeEntry(os, "relax", relax_);
170 }
171 
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 namespace Foam
176 {
178  (
181  );
182 }
183 
184 
185 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Run-time selectable general function of one variable.
Definition: Function1.H:64
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual const volScalarField & T() const =0
Temperature [K].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual tmp< scalarField > kappaEff(const label patchi) const =0
Effective thermal turbulent conductivity.
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:238
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:41
Uniform fixed heat flux boundary condition for Eulerian multi-phase cases. Constructs a mixed constra...
uniformFixedMultiphaseHeatFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
labelList f(nPoints)
dictionary dict
volScalarField & p