fixedTemperatureConstraint.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-2023 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 
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  );
45  }
46 
47  template<>
49  names[] =
50  {
51  "uniform",
52  "lookup"
53  };
54 }
55 
58 
59 
60 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
61 
62 void Foam::fv::fixedTemperatureConstraint::readCoeffs()
63 {
64  mode_ = modeNames_.read(coeffs().lookup("mode"));
65 
66  switch (mode_)
67  {
69  {
70  TValue_.reset
71  (
72  Function1<scalar>::New("temperature", coeffs()).ptr()
73  );
74  break;
75  }
77  {
78  TName_ = coeffs().lookupOrDefault<word>("T", "T");
79  break;
80  }
81  }
82 
83  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
84 
85  fraction_ =
86  coeffs().found("fraction")
87  ? Function1<scalar>::New("fraction", coeffs())
88  : autoPtr<Function1<scalar>>();
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93 
95 (
96  const word& name,
97  const word& modelType,
98  const fvMesh& mesh,
99  const dictionary& dict
100 )
101 :
102  fvConstraint(name, modelType, mesh, dict),
103  set_(mesh, coeffs()),
104  mode_(temperatureMode::uniform),
105  TValue_(nullptr),
106  TName_(word::null),
107  phaseName_(word::null)
108 {
109  readCoeffs();
110 }
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 {
117  const basicThermo& thermo =
118  mesh().lookupObject<basicThermo>
119  (
120  IOobject::groupName(physicalProperties::typeName, phaseName_)
121  );
122 
123  return wordList(1, thermo.he().name());
124 }
125 
126 
128 (
129  fvMatrix<scalar>& eqn,
130  const word& fieldName
131 ) const
132 {
133  const labelUList cells = set_.cells();
134 
135  const basicThermo& thermo =
136  mesh().lookupObject<basicThermo>
137  (
138  IOobject::groupName(physicalProperties::typeName, phaseName_)
139  );
140 
141  const scalar t = mesh().time().userTimeValue();
142 
143  switch (mode_)
144  {
145  case temperatureMode::uniform:
146  {
147  const scalarField Tuni(cells.size(), TValue_->value(t));
148  const scalarField heuni(thermo.he(Tuni, cells));
149 
150  if (fraction_.valid())
151  {
152  eqn.setValues
153  (
154  cells,
155  heuni,
156  scalarList(cells.size(), fraction_->value(t))
157  );
158  }
159  else
160  {
161  eqn.setValues(cells, heuni);
162  }
163 
164  break;
165  }
166  case temperatureMode::lookup:
167  {
168  const volScalarField& T =
169  mesh().lookupObject<volScalarField>(TName_);
170  const scalarField Tlkp(T, cells);
171  const scalarField helkp(thermo.he(Tlkp, cells));
172 
173  if (fraction_.valid())
174  {
175  eqn.setValues
176  (
177  cells,
178  helkp,
179  scalarList(cells.size(), fraction_->value(t))
180  );
181  }
182  else
183  {
184  eqn.setValues(cells, helkp);
185  }
186 
187  break;
188  }
189  }
190 
191  return cells.size();
192 }
193 
194 
196 {
197  set_.movePoints();
198  return true;
199 }
200 
201 
203 (
204  const polyTopoChangeMap& map
205 )
206 {
207  set_.topoChange(map);
208 }
209 
210 
212 {
213  set_.mapMesh(map);
214 }
215 
216 
218 (
219  const polyDistributionMap& map
220 )
221 {
222  set_.distribute(map);
223 }
224 
225 
227 {
229  {
230  set_.read(coeffs());
231  readCoeffs();
232  return true;
233  }
234  else
235  {
236  return false;
237  }
238 }
239 
240 
241 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
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:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
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
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:101
Fixed temperature equation constraint.
virtual bool movePoints()
Update for mesh motion.
fixedTemperatureConstraint(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
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.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
static const NamedEnum< temperatureMode, 2 > modeNames_
String representation of mode enums.
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)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31