lumpedMassTemperatureFvPatchScalarField.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) 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 
27 #include "fieldMapper.H"
29 #include "ZeroConstant.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 bool Foam::lumpedMassTemperatureFvPatchScalarField::closed() const
35 {
36  return mag(gSum(patch().Sf()))/gSum(patch().magSf()) < rootSmall;
37 }
38 
39 Foam::scalar Foam::lumpedMassTemperatureFvPatchScalarField::V() const
40 {
41  return -gSum(patch().Sf() & patch().Cf())
42  /patch().boundaryMesh().mesh().nSolutionD();
43 }
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
55  fixedValueFvPatchScalarField(p, iF),
56  rho_(dict.lookup<scalar>("rho")),
57  Cv_(dict.lookup<scalar>("Cv")),
58  T_
59  (
60  IOobject
61  (
62  "T_" + patch().name(),
63  db().time().name(),
64  db()
65  ),
66  dimensionedScalar(dimTemperature, dict.lookup<scalar>("T"))
67  ),
68  Q_
69  (
70  dict.found("Q")
71  ? Function1<scalar>::New
72  (
73  "Q",
74  db().time().userUnits(),
75  dimPower,
76  dict
77  )
78  : autoPtr<Function1<scalar>>(new Function1s::ZeroConstant<scalar>("Q"))
79  ),
80  V_(NaN)
81 {
82  if (!dict.readIfPresent("volume", V_))
83  {
84  if (closed())
85  {
86  V_ = V();
87  Info<< "Volume for the thermal mass, enclosed by patch '"
88  << patch().name() << "', = " << V_;
89  }
90  else
91  {
93  << "Patch '" << patch().name()
94  << "' corresponding to a thermal mass is not closed." << nl
95  << "Please specify the volume with the optional 'volume' entry."
96  << exit(FatalError);
97  }
98  }
99 
101 }
102 
103 
106 (
108  const fvPatch& p,
110  const fieldMapper& mapper
111 )
112 :
113  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
114  rho_(ptf.rho_),
115  Cv_(ptf.Cv_),
116  T_(ptf.T_),
117  Q_(ptf.Q_),
118  V_(ptf.V_)
119 {}
120 
121 
124 (
127 )
128 :
129  fixedValueFvPatchScalarField(ptf, iF),
130  rho_(ptf.rho_),
131  Cv_(ptf.Cv_),
132  T_(ptf.T_),
133  Q_(ptf.Q_),
134  V_(ptf.V_)
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 (
142  const fvPatchScalarField& ptf,
143  const fieldMapper& mapper
144 )
145 {
146  fixedValueFvPatchScalarField::map(ptf, mapper);
147 }
148 
149 
151 (
152  const fvPatchScalarField& ptf
153 )
154 {
155  fixedValueFvPatchScalarField::reset(ptf);
156 }
157 
158 
160 {
161  if
162  (
163  updated()
164  )
165  {
166  return;
167  }
168 
169  const thermophysicalTransportModel& ttm =
170  db().lookupType<thermophysicalTransportModel>
171  (
172  internalField().group()
173  );
174 
175  const label patchi = patch().index();
176  const scalarField Hf
177  (
178  ttm.kappaEff(patchi)*patch().magSf()*patch().deltaCoeffs()
179  );
180  const scalar Hs = rho_*Cv_*V_/db().time().deltaTValue();
181 
182  const scalar t = db().time().userTimeValue();
183 
184  T_.value() =
185  (
186  Q_->value(t) + gSum(Hf*patchInternalField())
187  + Hs*T_.oldTime().value()
188  )/(Hs + gSum(Hf));
189 
190  operator==(T_.value());
191 
192  fixedValueFvPatchScalarField::updateCoeffs();
193 }
194 
195 
197 (
198  Ostream& os
199 ) const
200 {
202  writeEntry(os, "rho", rho_);
203  writeEntry(os, "Cv", Cv_);
204  writeEntry(os, "T", T_.value());
205  writeEntry(os, Q_());
206  writeEntry(os, "volume", V_);
207  writeEntry(os, "value", *this);
208 }
209 
210 
211 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
212 
213 namespace Foam
214 {
216  (
219  );
220 }
221 
222 // ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Run-time selectable general function of one variable.
Definition: Function1.H:125
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Type & value()
Return a reference to the value.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition is applied to a patch which bounds a solid body, wholly or partially....
lumpedMassTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Abstract base class for all fluid and solid thermophysical transport models.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimPower
messageStream Info
const dimensionSet dimTemperature
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:266
dictionary dict
volScalarField & p