Loading...
Searching...
No Matches
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-2022 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"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace fv
39{
42}
43}
44
45
46// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47
49{
52 writeCommented(os, "Time");
53 writeTabbed(os, "nDampedCellsMin_[count]");
54 writeTabbed(os, "nDampedCellsMin_[%]");
55 writeTabbed(os, "nDampedCellsMax_[count]");
56 writeTabbed(os, "nDampedCellsMax_[%]");
57
58 os << endl;
60 writtenHeader_ = true;
61}
62
63
64// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65
67(
68 const word& name,
69 const word& modelType,
70 const dictionary& dict,
71 const fvMesh& mesh
72)
73:
74 fv::cellSetOption(name, modelType, dict, mesh),
75 writeFile(mesh, name, typeName, dict, false),
76 Tmin_(0),
77 Tmax_(0),
78 phase_()
79{
80 if (isActive())
81 {
83 }
84}
85
86
87// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88
90{
91 if (!(fv::cellSetOption::read(dict) && writeFile::read(dict)))
92 {
93 return false;
94 }
95
96 coeffs_.readEntry("min", Tmin_);
97 coeffs_.readEntry("max", Tmax_);
98 coeffs_.readIfPresent("phase", phase_);
99
100 if (Tmax_ < Tmin_)
101 {
103 << "Minimum temperature limit cannot exceed maximum limit" << nl
104 << "min = " << Tmin_ << nl
105 << "max = " << Tmax_
106 << exit(FatalIOError);
107 }
108
109 if (Tmin_ < 0)
110 {
112 << "Minimum temperature limit cannot be negative" << nl
113 << "min = " << Tmin_
114 << exit(FatalIOError);
115 }
116
117 // Set the field name to that of the energy
118 // field from which the temperature is obtained
119 const auto& thermo =
120 mesh_.lookupObject<basicThermo>
121 (
123 );
124
125 fieldNames_.resize(1, thermo.he().name());
126
128
129
130 if (canResetFile())
131 {
132 resetFile(typeName);
133 }
134
135 if (canWriteHeader())
136 {
137 writeFileHeader(file());
139
140
141 return true;
142}
143
144
146{
147 const auto& thermo =
148 mesh_.lookupObject<basicThermo>
149 (
151 );
152
153 scalarField Tmin(cells_.size(), Tmin_);
154 scalarField Tmax(cells_.size(), Tmax_);
155
156 scalarField heMin(thermo.he(thermo.p(), Tmin, cells_));
157 scalarField heMax(thermo.he(thermo.p(), Tmax, cells_));
158
159 scalarField& hec = he.primitiveFieldRef();
160
161 const scalarField& T = thermo.T();
162
163 scalar Tmin0 = min(T);
164 scalar Tmax0 = max(T);
165
166 // Count nTotCells ourselves
167 // (maybe only applying on a subset)
168 label nBelowMin(0);
169 label nAboveMax(0);
170 const label nTotCells(returnReduce(cells_.size(), sumOp<label>()));
171
172 forAll(cells_, i)
173 {
174 const label celli = cells_[i];
175 if (hec[celli] < heMin[i])
176 {
177 hec[celli] = heMin[i];
178 ++nBelowMin;
179 }
180 else if (hec[celli] > heMax[i])
181 {
182 hec[celli] = heMax[i];
183 ++nAboveMax;
184 }
185 }
186
187 reduce(nBelowMin, sumOp<label>());
188 reduce(nAboveMax, sumOp<label>());
189
190 reduce(Tmin0, minOp<scalar>());
191 reduce(Tmax0, maxOp<scalar>());
192
193 // Percent, max 2 decimal places
194 const auto percent = [](scalar num, label denom) -> scalar
195 {
196 return (denom ? 1e-2*round(1e4*num/denom) : 0);
197 };
198
199 const scalar nBelowMinPercent = percent(nBelowMin, nTotCells);
200 const scalar nAboveMaxPercent = percent(nAboveMax, nTotCells);
201
202 Info<< type() << "=" << name_ << ", Type=Lower"
203 << ", LimitedCells=" << nBelowMin
204 << ", CellsPercent=" << nBelowMinPercent
205 << ", Tmin=" << Tmin_
206 << ", UnlimitedTmin=" << Tmin0
207 << endl;
208
209 Info<< type() << "=" << name_ << ", Type=Upper"
210 << ", LimitedCells=" << nAboveMax
211 << ", CellsPercent=" << nAboveMaxPercent
212 << ", Tmax=" << Tmax_
213 << ", UnlimitedTmax=" << Tmax0
214 << endl;
215
216
217 if (canWriteToFile())
218 {
219 file()
220 << mesh_.time().timeOutputValue() << token::TAB
221 << nBelowMin << token::TAB
222 << nBelowMinPercent << token::TAB
223 << nAboveMax << token::TAB
224 << nAboveMaxPercent
225 << endl;
226 }
227
228
229 // Handle boundaries in the case of 'all'
230 bool changedValues = (nBelowMin || nAboveMax);
232 {
233 volScalarField::Boundary& bf = he.boundaryFieldRef();
234
235 forAll(bf, patchi)
236 {
237 fvPatchScalarField& hep = bf[patchi];
238
239 if (!hep.fixesValue())
240 {
241 const scalarField& pp = thermo.p().boundaryField()[patchi];
242
243 scalarField Tminp(pp.size(), Tmin_);
244 scalarField Tmaxp(pp.size(), Tmax_);
245
246 scalarField heMinp(thermo.he(pp, Tminp, patchi));
247 scalarField heMaxp(thermo.he(pp, Tmaxp, patchi));
248
249 forAll(hep, facei)
250 {
251 if (hep[facei] < heMinp[facei])
252 {
253 hep[facei] = heMinp[facei];
254 changedValues = true;
255 }
256 else if (hep[facei] > heMaxp[facei])
257 {
258 hep[facei] = heMaxp[facei];
259 changedValues = true;
260 }
261 }
262 }
263 }
264 }
265
266
267 if (returnReduceOr(changedValues))
268 {
269 // We've changed internal values so give
270 // boundary conditions opportunity to correct
271 he.correctBoundaryConditions();
272 }
273}
274
275
276// ************************************************************************* //
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.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
virtual bool canWriteToFile() const
Flag to allow writing to the file.
Definition writeFile.C:292
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
bool writtenHeader_
Flag to identify whether the header has been written.
Definition writeFile.H:157
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void resetFile(const word &name)
Reset internal file pointer to new file with new name.
Definition writeFile.C:165
virtual bool canWriteHeader() const
Flag to allow writing the header.
Definition writeFile.C:304
virtual bool canResetFile() const
Flag to allow resetting the file.
Definition writeFile.C:298
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual bool fixesValue() const
True if the patch field fixes a value.
Intermediate abstract class for handling cell-set options for the derived fvOptions.
labelList cells_
Set of cells to apply source to.
virtual bool read(const dictionary &dict)
Read source dictionary.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
bool useSubMesh() const noexcept
True if sub-selection should be used.
virtual bool isActive()
Is the source active?
Corrects temperature field (i.e. T) within a specified region by applying limits between a given mini...
word phase_
Optional phase name.
virtual bool read(const dictionary &dict)
Read dictionary.
limitTemperature(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
scalar Tmax_
Maximum temperature limit [K].
virtual void correct(volScalarField &he)
Correct the energy field.
scalar Tmin_
Minimum temperature limit [K].
void writeFileHeader(Ostream &os)
Write file header information.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
@ TAB
Tab [isspace].
Definition token.H:142
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
Namespace for finite-volume.
Namespace for OpenFOAM.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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 POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
psiReactionThermo & thermo
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299