limitTemperature.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "limitTemperature.H"
30#include "fvMesh.H"
31#include "basicThermo.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace fv
39{
42}
43}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
49(
50 const word& name,
51 const word& modelType,
52 const dictionary& dict,
53 const fvMesh& mesh
54)
55:
56 fv::cellSetOption(name, modelType, dict, mesh),
57 Tmin_(coeffs_.get<scalar>("min")),
58 Tmax_(coeffs_.get<scalar>("max")),
59 phase_(coeffs_.getOrDefault<word>("phase", word::null))
60{
61 // Set the field name to that of the energy
62 // field from which the temperature is obtained
63 const auto& thermo =
65 (
67 );
68
69 fieldNames_.resize(1, thermo.he().name());
70
72}
73
74
75// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76
78{
80 {
81 coeffs_.readEntry("min", Tmin_);
82 coeffs_.readEntry("max", Tmax_);
83
84 if (Tmax_ < Tmin_)
85 {
87 << "Minimum temperature limit cannot exceed maximum limit" << nl
88 << "min = " << Tmin_ << nl
89 << "max = " << Tmax_
91 }
92
93 if (Tmin_ < 0)
94 {
96 << "Minimum temperature limit cannot be negative" << nl
97 << "min = " << Tmin_
99 }
100
101 return true;
102 }
103
104 return false;
105}
106
107
109{
110 const auto& thermo =
111 mesh_.lookupObject<basicThermo>
112 (
114 );
115
116 scalarField Tmin(cells_.size(), Tmin_);
117 scalarField Tmax(cells_.size(), Tmax_);
118
119 scalarField heMin(thermo.he(thermo.p(), Tmin, cells_));
120 scalarField heMax(thermo.he(thermo.p(), Tmax, cells_));
121
122 scalarField& hec = he.primitiveFieldRef();
123
124 const scalarField& T = thermo.T();
125
126 scalar Tmin0 = min(T);
127 scalar Tmax0 = max(T);
128
129 // Count nTotCells ourselves
130 // (maybe only applying on a subset)
131 label nBelowMin(0);
132 label nAboveMax(0);
133 const label nTotCells(returnReduce(cells_.size(), sumOp<label>()));
134
135 forAll(cells_, i)
136 {
137 const label celli = cells_[i];
138 if (hec[celli] < heMin[i])
139 {
140 hec[celli] = heMin[i];
141 ++nBelowMin;
142 }
143 else if (hec[celli] > heMax[i])
144 {
145 hec[celli] = heMax[i];
146 ++nAboveMax;
147 }
148 }
149
150 reduce(nBelowMin, sumOp<label>());
151 reduce(nAboveMax, sumOp<label>());
152
153 reduce(Tmin0, minOp<scalar>());
154 reduce(Tmax0, maxOp<scalar>());
155
156 // Percent, max 2 decimal places
157 const auto percent = [](scalar num, label denom) -> scalar
158 {
159 return (denom ? 1e-2*round(1e4*num/denom) : 0);
160 };
161
162 Info<< type() << ' ' << name_ << " Lower limited " << nBelowMin << " ("
163 << percent(nBelowMin, nTotCells)
164 << "%) of cells, with min limit " << Tmin_ << endl;
165
166 Info<< type() << ' ' << name_ << " Upper limited " << nAboveMax << " ("
167 << percent(nAboveMax, nTotCells)
168 << "%) of cells, with max limit " << Tmax_ << endl;
169
170 Info<< type() << ' ' << name_ << " Unlimited Tmin " << Tmin0 << endl;
171 Info<< type() << ' ' << name_ << " Unlimited Tmax " << Tmax0 << endl;
172
173
174 // Handle boundaries in the case of 'all'
175 bool changedValues = (nBelowMin || nAboveMax);
177 {
178 volScalarField::Boundary& bf = he.boundaryFieldRef();
179
180 forAll(bf, patchi)
181 {
182 fvPatchScalarField& hep = bf[patchi];
183
184 if (!hep.fixesValue())
185 {
186 const scalarField& pp = thermo.p().boundaryField()[patchi];
187
188 scalarField Tminp(pp.size(), Tmin_);
189 scalarField Tmaxp(pp.size(), Tmax_);
190
191 scalarField heMinp(thermo.he(pp, Tminp, patchi));
192 scalarField heMaxp(thermo.he(pp, Tmaxp, patchi));
193
194 forAll(hep, facei)
195 {
196 if (hep[facei] < heMinp[facei])
197 {
198 hep[facei] = heMinp[facei];
199 changedValues = true;
200 }
201 else if (hep[facei] > heMaxp[facei])
202 {
203 hep[facei] = heMaxp[facei];
204 changedValues = true;
205 }
206 }
207 }
208 }
209 }
210
211
212 if (returnReduce(changedValues, orOp<bool>()))
213 {
214 // We've changed internal values so give
215 // boundary conditions opportunity to correct
216 he.correctBoundaryConditions();
217 }
218}
219
220
221// ************************************************************************* //
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.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read()
Re-read model coefficients if they have changed.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:66
static const word dictName
Definition: basicThermo.H:256
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:337
Intermediate abstract class for handling cell-set options for the derived fvOptions.
bool useSubMesh() const noexcept
True if sub-selection should be used.
Corrects temperature field (i.e. T) within a specified region by applying limits between a given mini...
word phase_
Optional phase name [K].
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:127
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:148
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:48
const Type & lookupObject(const word &name, const bool recursive=false) const
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
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
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
IOerror FatalIOError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
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
labelList fv(nPoints)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333