fixedTemperature.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-2025 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 "fixedTemperature.H"
27 #include "fvMesh.H"
28 #include "fvMatrices.H"
29 #include "basicThermo.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
40  (
44  );
46  (
49  dictionary,
50  fixedTemperatureConstraint,
51  "fixedTemperatureConstraint"
52  );
53 }
54 }
55 
56 
59 {
60  "uniform",
61  "lookup"
62 };
63 
64 
65 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
66 
67 void Foam::fv::fixedTemperature::readCoeffs(const dictionary& dict)
68 {
69  mode_ = modeNames_.read(dict.lookup("mode"));
70 
71  switch (mode_)
72  {
74  {
75  TValue_.reset
76  (
78  (
79  "temperature",
80  mesh().time().userUnits(),
82  dict
83  ).ptr()
84  );
85  break;
86  }
88  {
89  TName_ = dict.lookupOrDefault<word>("T", "T");
90  break;
91  }
92  }
93 
94  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
95 
96  fraction_ =
97  dict.found("fraction")
99  (
100  "fraction",
101  mesh().time().userUnits(),
102  unitFraction,
103  dict
104  )
105  : autoPtr<Function1<scalar>>();
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110 
112 (
113  const word& name,
114  const word& modelType,
115  const fvMesh& mesh,
116  const dictionary& dict
117 )
118 :
119  fvConstraint(name, modelType, mesh, dict),
120  zone_(mesh, coeffs(dict)),
121  mode_(temperatureMode::uniform),
122  TValue_(nullptr),
123  TName_(word::null),
124  phaseName_(word::null)
125 {
126  readCoeffs(coeffs(dict));
127 }
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
133 {
134  const basicThermo& thermo =
136  (
137  IOobject::groupName(physicalProperties::typeName, phaseName_)
138  );
139 
140  return wordList(1, thermo.he().name());
141 }
142 
143 
145 (
146  fvMatrix<scalar>& eqn,
147  const word& fieldName
148 ) const
149 {
150  const labelList& cells = zone_.zone();
151 
152  const basicThermo& thermo =
154  (
155  IOobject::groupName(physicalProperties::typeName, phaseName_)
156  );
157 
158  const scalar t = mesh().time().value();
159 
160  switch (mode_)
161  {
162  case temperatureMode::uniform:
163  {
164  const scalarField Tuni(cells.size(), TValue_->value(t));
165  const scalarField heuni(thermo.he(Tuni, cells));
166 
167  if (fraction_.valid())
168  {
169  eqn.setValues
170  (
171  cells,
172  heuni,
173  scalarList(cells.size(), fraction_->value(t))
174  );
175  }
176  else
177  {
178  eqn.setValues(cells, heuni);
179  }
180 
181  break;
182  }
183  case temperatureMode::lookup:
184  {
185  const volScalarField& T =
186  mesh().lookupObject<volScalarField>(TName_);
187  const scalarField Tlkp(T, cells);
188  const scalarField helkp(thermo.he(Tlkp, cells));
189 
190  if (fraction_.valid())
191  {
192  eqn.setValues
193  (
194  cells,
195  helkp,
196  scalarList(cells.size(), fraction_->value(t))
197  );
198  }
199  else
200  {
201  eqn.setValues(cells, helkp);
202  }
203 
204  break;
205  }
206  }
207 
208  return cells.size();
209 }
210 
211 
213 {
214  zone_.movePoints();
215  return true;
216 }
217 
218 
220 {
221  zone_.topoChange(map);
222 }
223 
224 
226 {
227  zone_.mapMesh(map);
228 }
229 
230 
232 {
233  zone_.distribute(map);
234 }
235 
236 
238 {
240  {
241  zone_.read(coeffs(dict));
242  readCoeffs(coeffs(dict));
243  return true;
244  }
245  else
246  {
247  return false;
248  }
249 }
250 
251 
252 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< scalar > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricField class.
static word groupName(Name name, const word &group)
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:54
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
Finite volume options abstract base class.
Definition: fvConstraint.H:57
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvConstraint.C:164
const dictionary & coeffs(const dictionary &) const
Return the coefficients sub-dictionary.
Definition: fvConstraintI.H:41
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvConstraintI.H:34
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void setValues(const labelUList &cells, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:502
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
Fixed temperature equation constraint.
virtual bool movePoints()
Update for mesh motion.
static const NamedEnum< temperatureMode, 2 > modeNames_
String representation of mode enums.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual bool constrain(fvMatrix< scalar > &eqn, const word &fieldName) const
Constrain energy equation to fix the temperature.
fixedTemperature(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
temperatureMode
Temperature mode.
virtual wordList constrainedFields() const
Return the list of fields constrained by the fvConstraint.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const cellShapeList & cells
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
addBackwardCompatibleToRunTimeSelectionTable(fvConstraint, fixedTemperature, dictionary, fixedTemperatureConstraint, "fixedTemperatureConstraint")
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTemperature
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const unitConversion unitFraction
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31