limitTemperature.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) 2012-2022 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 
26 #include "limitTemperature.H"
27 #include "basicThermo.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
36  defineTypeNameAndDebug(limitTemperature, 0);
38  (
39  fvConstraint,
40  limitTemperature,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::fv::limitTemperature::readCoeffs()
50 {
51  Tmin_ = coeffs().lookup<scalar>("min");
52  Tmax_ = coeffs().lookup<scalar>("max");
53  fieldName_ = coeffs().lookupOrDefault<word>("field", word::null);
54  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const word& name,
63  const word& modelType,
64  const dictionary& dict,
65  const fvMesh& mesh
66 )
67 :
68  fvConstraint(name, modelType, dict, mesh),
69  set_(coeffs(), mesh),
70  Tmin_(-vGreat),
71  Tmax_(vGreat),
72  fieldName_(word::null),
73  phaseName_(word::null)
74 {
75  readCoeffs();
76 }
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
82 {
83  if (fieldName_ != word::null)
84  {
85  return wordList(1, IOobject::groupName(fieldName_, phaseName_));
86  }
87  else
88  {
89  const basicThermo& thermo =
90  mesh().lookupObject<basicThermo>
91  (
92  IOobject::groupName(physicalProperties::typeName, phaseName_)
93  );
94 
95  return wordList(1, thermo.he().name());
96  }
97 }
98 
99 
101 {
102  const labelList& cells = set_.cells();
103 
104  if (he.dimensions() == dimTemperature)
105  {
106  scalarField& Tc = he.primitiveFieldRef();
107 
108  forAll(cells, i)
109  {
110  const label celli = cells[i];
111  Tc[celli] = max(min(Tc[celli], Tmax_), Tmin_);
112  }
113 
114  // Handle boundaries in the case of 'all'
115  if (set_.selectionMode() == fvCellSet::selectionModeType::all)
116  {
118 
119  forAll(Tbf, patchi)
120  {
121  fvPatchScalarField& Tp = Tbf[patchi];
122 
123  if (!Tp.fixesValue())
124  {
125  forAll(Tp, facei)
126  {
127  Tp[facei] = max(min(Tp[facei], Tmax_), Tmin_);
128  }
129  }
130  }
131  }
132  }
133  else
134  {
135  const basicThermo& thermo =
136  mesh().lookupObject<basicThermo>
137  (
138  IOobject::groupName(physicalProperties::typeName, phaseName_)
139  );
140 
141  const scalarField Tmin(cells.size(), Tmin_);
142  const scalarField Tmax(cells.size(), Tmax_);
143 
144  const scalarField heMin(thermo.he(Tmin, cells));
145  const scalarField heMax(thermo.he(Tmax, cells));
146 
147  scalarField& hec = he.primitiveFieldRef();
148 
149  forAll(cells, i)
150  {
151  const label celli = cells[i];
152  hec[celli] = max(min(hec[celli], heMax[i]), heMin[i]);
153  }
154 
155  // Handle boundaries in the case of 'all'
156  if (set_.selectionMode() == fvCellSet::selectionModeType::all)
157  {
159 
160  forAll(bf, patchi)
161  {
162  fvPatchScalarField& hep = bf[patchi];
163 
164  if (!hep.fixesValue())
165  {
166  const scalarField Tminp(hep.size(), Tmin_);
167  const scalarField Tmaxp(hep.size(), Tmax_);
168 
169  const scalarField heMinp(thermo.he(Tminp, patchi));
170  const scalarField heMaxp(thermo.he(Tmaxp, patchi));
171 
172  forAll(hep, facei)
173  {
174  hep[facei] =
175  max(min(hep[facei], heMaxp[facei]), heMinp[facei]);
176  }
177  }
178  }
179  }
180  }
181 
182  return cells.size();
183 }
184 
185 
187 {
188  set_.movePoints();
189  return true;
190 }
191 
192 
194 {
195  set_.topoChange(map);
196 }
197 
198 
200 {
201  set_.mapMesh(map);
202 }
203 
204 
206 {
207  set_.distribute(map);
208 }
209 
210 
212 {
213  if (fvConstraint::read(dict))
214  {
215  set_.read(coeffs());
216  readCoeffs();
217  return true;
218  }
219  else
220  {
221  return false;
222  }
223 }
224 
225 
226 // ************************************************************************* //
limitTemperature(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const word & name() const
Return name.
Definition: IOobject.H:315
fluidReactionThermo & thermo
Definition: createFields.H:28
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:77
virtual bool movePoints()
Update for mesh motion.
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:310
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual bool constrain(volScalarField &he) const
Constrain the energy field.
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
const dimensionSet & dimensions() const
Return dimensions.
const cellShapeList & cells
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static word groupName(Name name, const word &group)
static const word null
An empty word.
Definition: word.H:77
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
thermo he()
virtual wordList constrainedFields() const
Return the list of fields constrained by the fvConstraint.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvConstraint.C:143
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
Finite volume options abstract base class.
Definition: fvConstraint.H:56
const dimensionSet dimTemperature
Namespace for OpenFOAM.