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 -------------------------------------------------------------------------------
11 License
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 
36 namespace Foam
37 {
38 namespace fv
39 {
40  defineTypeNameAndDebug(limitTemperature, 0);
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;
59 
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 {
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());
138  }
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_ << " Lower limited " << nBelowMin << " ("
203  << nBelowMinPercent
204  << "%) of cells, with min limit " << Tmin_ << endl;
205 
206  Info<< type() << ' ' << name_ << " Upper limited " << nAboveMax << " ("
207  << nAboveMaxPercent
208  << "%) of cells, with max limit " << Tmax_ << endl;
209 
210  Info<< type() << ' ' << name_ << " Unlimited Tmin " << Tmin0 << endl;
211  Info<< type() << ' ' << name_ << " Unlimited Tmax " << Tmax0 << endl;
212 
213 
214  if (canWriteToFile())
215  {
216  file()
217  << mesh_.time().timeOutputValue() << token::TAB
218  << nBelowMin << token::TAB
219  << nBelowMinPercent << token::TAB
220  << nAboveMax << token::TAB
221  << nAboveMaxPercent
222  << endl;
223  }
224 
225 
226  // Handle boundaries in the case of 'all'
227  bool changedValues = (nBelowMin || nAboveMax);
229  {
230  volScalarField::Boundary& bf = he.boundaryFieldRef();
231 
232  forAll(bf, patchi)
233  {
234  fvPatchScalarField& hep = bf[patchi];
235 
236  if (!hep.fixesValue())
237  {
238  const scalarField& pp = thermo.p().boundaryField()[patchi];
239 
240  scalarField Tminp(pp.size(), Tmin_);
241  scalarField Tmaxp(pp.size(), Tmax_);
242 
243  scalarField heMinp(thermo.he(pp, Tminp, patchi));
244  scalarField heMaxp(thermo.he(pp, Tmaxp, patchi));
245 
246  forAll(hep, facei)
247  {
248  if (hep[facei] < heMinp[facei])
249  {
250  hep[facei] = heMinp[facei];
251  changedValues = true;
252  }
253  else if (hep[facei] > heMaxp[facei])
254  {
255  hep[facei] = heMaxp[facei];
256  changedValues = true;
257  }
258  }
259  }
260  }
261  }
262 
263 
264  if (returnReduceOr(changedValues))
265  {
266  // We've changed internal values so give
267  // boundary conditions opportunity to correct
268  he.correctBoundaryConditions();
269  }
270 }
271 
272 
273 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
limitTemperature(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
dictionary dict
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Corrects temperature field (i.e. T) within a specified region by applying limits between a given mini...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Tab [isspace].
Definition: token.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual void correct(volScalarField &he)
Correct the energy field.
virtual bool read(const dictionary &dict)
Read source dictionary.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
bool writtenHeader_
Flag to identify whether the header has been written.
Definition: writeFile.H:157
void writeFileHeader(Ostream &os)
Write file header information.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
psiReactionThermo & thermo
Definition: createFields.H:28
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
scalar Tmin_
Minimum temperature limit [K].
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual bool isActive()
Is the source active?
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
scalar Tmax_
Maximum temperature limit [K].
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
OBJstream os(runTime.globalPath()/outputName)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Intermediate abstract class for handling cell-set options for the derived fvOptions.
virtual bool read(const dictionary &dict)
Read dictionary.
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
bool useSubMesh() const noexcept
True if sub-selection should be used.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:123
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329