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-2020 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  {
38  defineTypeNameAndDebug(fixedTemperatureConstraint, 0);
40  (
41  option,
42  fixedTemperatureConstraint,
43  dictionary
44  );
45  }
46 
47  template<>
48  const char* NamedEnum<fv::fixedTemperatureConstraint::temperatureMode, 2>::
49  names[] =
50  {
51  "uniform",
52  "lookup"
53  };
54 }
55 
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const word& name,
65  const word& modelType,
66  const dictionary& dict,
67  const fvMesh& mesh
68 )
69 :
70  cellSetOption(name, modelType, dict, mesh),
71  mode_(temperatureModeNames_.read(coeffs_.lookup("mode"))),
72  Tuniform_(nullptr),
73  TName_("T"),
74  phase_(coeffs_.lookupOrDefault<word>("phase", word::null))
75 {
76  switch (mode_)
77  {
78  case tmUniform:
79  {
80  Tuniform_.reset
81  (
82  Function1<scalar>::New("temperature", coeffs_).ptr()
83  );
84  break;
85  }
86  case tmLookup:
87  {
88  TName_ = coeffs_.lookupOrDefault<word>("T", "T");
89  break;
90  }
91  default:
92  {
93  // error handling done by NamedEnum
94  }
95  }
96 
97  // Set the field name to that of the energy field from which the temperature
98  // is obtained
99 
100  const basicThermo& thermo =
101  mesh_.lookupObject<basicThermo>
102  (
104  );
105 
106  fieldNames_.setSize(1, thermo.he().name());
107 
108  applied_.setSize(1, false);
109 }
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 (
116  fvMatrix<scalar>& eqn,
117  const label
118 )
119 {
120  const basicThermo& thermo =
121  mesh_.lookupObject<basicThermo>
122  (
124  );
125 
126  switch (mode_)
127  {
128  case tmUniform:
129  {
130  const scalar t = mesh_.time().value();
131  scalarField Tuni(cells_.size(), Tuniform_->value(t));
132  eqn.setValues(cells_, thermo.he(Tuni, cells_));
133 
134  break;
135  }
136  case tmLookup:
137  {
138  const volScalarField& T =
139  mesh().lookupObject<volScalarField>(TName_);
140 
141  scalarField Tlkp(T, cells_);
142  eqn.setValues(cells_, thermo.he(Tlkp, cells_));
143 
144  break;
145  }
146  default:
147  {
148  // error handling done by NamedEnum
149  }
150  }
151 }
152 
153 
155 {
156  if (cellSetOption::read(dict))
157  {
158  if (coeffs_.found(Tuniform_->name()))
159  {
160  Tuniform_.reset
161  (
162  Function1<scalar>::New(Tuniform_->name(), dict).ptr()
163  );
164  }
165 
166  coeffs_.readIfPresent("T", TName_);
167 
168  return true;
169  }
170  else
171  {
172  return false;
173  }
174 }
175 
176 
177 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:52
dictionary dict
virtual bool read(const dictionary &dict)
Read dictionary.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:303
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
fixedTemperatureConstraint(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
virtual bool read(const dictionary &dict)
Read source dictionary.
rhoReactionThermo & thermo
Definition: createFields.H:28
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Macros for easy insertion into run-time selection tables.
static const NamedEnum< temperatureMode, 2 > temperatureModeNames_
String representation of temperatureMode enums.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
dynamicFvMesh & mesh
virtual void constrain(fvMatrix< scalar > &eqn, const label fieldi)
Constrain energy equation to fix the temperature.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static word groupName(Name name, const word &group)
void setValues(const labelUList &cells, const UList< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:484
const Type & value() const
Return const reference to value.
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. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const Time & time() const
Return time.
Definition: IOobject.C:328
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Cell-set options abstract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:69
Namespace for OpenFOAM.