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-2024 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 
57 namespace Foam
58 {
59  template<>
61  {"uniform", "lookup"};
62 }
63 
66 
67 
68 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
69 
70 void Foam::fv::fixedTemperature::readCoeffs()
71 {
72  mode_ = modeNames_.read(coeffs().lookup("mode"));
73 
74  switch (mode_)
75  {
77  {
78  TValue_.reset
79  (
81  (
82  "temperature",
83  mesh().time().userUnits(),
85  coeffs()
86  ).ptr()
87  );
88  break;
89  }
91  {
92  TName_ = coeffs().lookupOrDefault<word>("T", "T");
93  break;
94  }
95  }
96 
97  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
98 
99  fraction_ =
100  coeffs().found("fraction")
102  (
103  "fraction",
104  mesh().time().userUnits(),
105  unitFraction,
106  coeffs()
107  )
108  : autoPtr<Function1<scalar>>();
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 
115 (
116  const word& name,
117  const word& modelType,
118  const fvMesh& mesh,
119  const dictionary& dict
120 )
121 :
122  fvConstraint(name, modelType, mesh, dict),
123  set_(mesh, coeffs()),
124  mode_(temperatureMode::uniform),
125  TValue_(nullptr),
126  TName_(word::null),
127  phaseName_(word::null)
128 {
129  readCoeffs();
130 }
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
136 {
137  const basicThermo& thermo =
138  mesh().lookupObject<basicThermo>
139  (
140  IOobject::groupName(physicalProperties::typeName, phaseName_)
141  );
142 
143  return wordList(1, thermo.he().name());
144 }
145 
146 
148 (
149  fvMatrix<scalar>& eqn,
150  const word& fieldName
151 ) const
152 {
153  const labelUList cells = set_.cells();
154 
155  const basicThermo& thermo =
156  mesh().lookupObject<basicThermo>
157  (
158  IOobject::groupName(physicalProperties::typeName, phaseName_)
159  );
160 
161  const scalar t = mesh().time().value();
162 
163  switch (mode_)
164  {
165  case temperatureMode::uniform:
166  {
167  const scalarField Tuni(cells.size(), TValue_->value(t));
168  const scalarField heuni(thermo.he(Tuni, cells));
169 
170  if (fraction_.valid())
171  {
172  eqn.setValues
173  (
174  cells,
175  heuni,
176  scalarList(cells.size(), fraction_->value(t))
177  );
178  }
179  else
180  {
181  eqn.setValues(cells, heuni);
182  }
183 
184  break;
185  }
186  case temperatureMode::lookup:
187  {
188  const volScalarField& T =
189  mesh().lookupObject<volScalarField>(TName_);
190  const scalarField Tlkp(T, cells);
191  const scalarField helkp(thermo.he(Tlkp, cells));
192 
193  if (fraction_.valid())
194  {
195  eqn.setValues
196  (
197  cells,
198  helkp,
199  scalarList(cells.size(), fraction_->value(t))
200  );
201  }
202  else
203  {
204  eqn.setValues(cells, helkp);
205  }
206 
207  break;
208  }
209  }
210 
211  return cells.size();
212 }
213 
214 
216 {
217  set_.movePoints();
218  return true;
219 }
220 
221 
223 {
224  set_.topoChange(map);
225 }
226 
227 
229 {
230  set_.mapMesh(map);
231 }
232 
233 
235 {
236  set_.distribute(map);
237 }
238 
239 
241 {
243  {
244  set_.read(coeffs());
245  readCoeffs();
246  return true;
247  }
248  else
249  {
250  return false;
251  }
252 }
253 
254 
255 // ************************************************************************* //
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:54
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:61
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
Finite volume options abstract base class.
Definition: fvConstraint.H:57
const dictionary & coeffs() const
Return dictionary.
Definition: fvConstraintI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvConstraint.C:166
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:99
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.
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
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimTemperature
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const unitConversion unitFraction
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31