limitTemperature.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2012-2016 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 
50 Foam::fv::limitTemperature::limitTemperature
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_(readScalar(coeffs_.lookup("min"))),
60  Tmax_(readScalar(coeffs_.lookup("max")))
61 {
62  // Set the field name to that of the energy field from which the temperature
63  // is obtained
64 
65  const basicThermo& thermo =
66  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
67 
68  fieldNames_.setSize(1, thermo.he().name());
69 
70  applied_.setSize(1, false);
71 }
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
78  if (cellSetOption::read(dict))
79  {
80  coeffs_.lookup("min") >> Tmin_;
81  coeffs_.lookup("max") >> Tmax_;
82 
83  return true;
84  }
85  else
86  {
87  return false;
88  }
89 }
90 
91 
93 {
94  const basicThermo& thermo =
95  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
96 
97  scalarField Tmin(cells_.size(), Tmin_);
98  scalarField Tmax(cells_.size(), Tmax_);
99 
100  scalarField heMin(thermo.he(thermo.p(), Tmin, cells_));
101  scalarField heMax(thermo.he(thermo.p(), Tmax, cells_));
102 
103  scalarField& hec = he.primitiveFieldRef();
104 
105  forAll(cells_, i)
106  {
107  label celli = cells_[i];
108  hec[celli]= max(min(hec[celli], heMax[i]), heMin[i]);
109  }
110 
111  // handle boundaries in the case of 'all'
112  if (selectionMode_ == smAll)
113  {
114  volScalarField::Boundary& bf = he.boundaryFieldRef();
115 
116  forAll(bf, patchi)
117  {
118  fvPatchScalarField& hep = bf[patchi];
119 
120  if (!hep.fixesValue())
121  {
122  const scalarField& pp = thermo.p().boundaryField()[patchi];
123 
124  scalarField Tminp(pp.size(), Tmin_);
125  scalarField Tmaxp(pp.size(), Tmax_);
126 
127  scalarField heMinp(thermo.he(pp, Tminp, patchi));
128  scalarField heMaxp(thermo.he(pp, Tmaxp, patchi));
129 
130  forAll(hep, facei)
131  {
132  hep[facei] =
133  max(min(hep[facei], heMaxp[facei]), heMinp[facei]);
134  }
135  }
136  }
137  }
138 }
139 
140 
141 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
volScalarField & he
Definition: YEEqn.H:51
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:291
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:315
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
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.
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:115
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:477
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
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 abtract 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.