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-2019 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 "fvMesh.H"
28 #include "basicThermo.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
37  defineTypeNameAndDebug(limitTemperature, 0);
39  (
40  option,
41  limitTemperature,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const word& name,
53  const word& modelType,
54  const dictionary& dict,
55  const fvMesh& mesh
56 )
57 :
58  cellSetOption(name, modelType, dict, mesh),
59  Tmin_(coeffs_.lookup<scalar>("min")),
60  Tmax_(coeffs_.lookup<scalar>("max")),
61  phase_(coeffs_.lookupOrDefault<word>("phase", word::null))
62 {
63  // Set the field name to that of the energy field from which the temperature
64  // is obtained
65  const basicThermo& thermo =
66  mesh_.lookupObject<basicThermo>
67  (
69  );
70 
71  fieldNames_.setSize(1, thermo.he().name());
72 
73  applied_.setSize(1, false);
74 }
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  if (cellSetOption::read(dict))
82  {
83  coeffs_.lookup("min") >> Tmin_;
84  coeffs_.lookup("max") >> Tmax_;
85 
86  return true;
87  }
88  else
89  {
90  return false;
91  }
92 }
93 
94 
96 {
97  const basicThermo& thermo =
98  mesh_.lookupObject<basicThermo>
99  (
101  );
102 
103  scalarField Tmin(cells_.size(), Tmin_);
104  scalarField Tmax(cells_.size(), Tmax_);
105 
106  scalarField heMin(thermo.he(Tmin, cells_));
107  scalarField heMax(thermo.he(Tmax, cells_));
108 
109  scalarField& hec = he.primitiveFieldRef();
110 
111  forAll(cells_, i)
112  {
113  label celli = cells_[i];
114  hec[celli]= max(min(hec[celli], heMax[i]), heMin[i]);
115  }
116 
117  // handle boundaries in the case of 'all'
118  if (selectionMode_ == smAll)
119  {
120  volScalarField::Boundary& bf = he.boundaryFieldRef();
121 
122  forAll(bf, patchi)
123  {
124  fvPatchScalarField& hep = bf[patchi];
125 
126  if (!hep.fixesValue())
127  {
128  scalarField Tminp(hep.size(), Tmin_);
129  scalarField Tmaxp(hep.size(), Tmax_);
130 
131  scalarField heMinp(thermo.he(Tminp, patchi));
132  scalarField heMaxp(thermo.he(Tmaxp, patchi));
133 
134  forAll(hep, facei)
135  {
136  hep[facei] =
137  max(min(hep[facei], heMaxp[facei]), heMinp[facei]);
138  }
139  }
140  }
141  }
142 }
143 
144 
145 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
volScalarField & he
Definition: YEEqn.H:50
limitTemperature(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:300
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual void correct(volScalarField &he)
Correct the energy field.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
virtual bool read(const dictionary &dict)
Read source dictionary.
rhoReactionThermo & thermo
Definition: createFields.H:28
Macros for easy insertion into run-time selection tables.
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
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static word groupName(Name name, const word &group)
static const word null
An empty word.
Definition: word.H:77
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Cell-set options abstract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:69
virtual bool read(const dictionary &dict)
Read dictionary.
Namespace for OpenFOAM.