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-2025 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 
40 Foam::scalar Foam::lumpedMassTemperatureFvPatchScalarField::V() const
41 {
42  return
43  -gSum(patch().Sf() & patch().Cf())
44  /patch().boundaryMesh().mesh().nSolutionD();
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
52 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
58  fixedValueFvPatchScalarField(p, iF),
59  rho_(dict.lookup<scalar>("rho")),
60  Cv_(dict.lookup<scalar>("Cv")),
61  T_
62  (
63  IOobject
64  (
65  "T_" + patch().name(),
66  db().time().name(),
67  db()
68  ),
69  dimensionedScalar(dimTemperature, dict.lookup<scalar>("T"))
70  ),
71  Q_
72  (
73  dict.found("Q")
74  ? Function1<scalar>::New
75  (
76  "Q",
77  db().time().userUnits(),
78  dimPower,
79  dict
80  )
81  : autoPtr<Function1<scalar>>(new Function1s::ZeroConstant<scalar>("Q"))
82  ),
83  V_(NaN)
84 {
85  if (!dict.readIfPresent("volume", V_))
86  {
87  if (closed())
88  {
89  V_ = V();
90  Info<< "Volume for the thermal mass, enclosed by patch '"
91  << patch().name() << "', = " << V_;
92  }
93  else
94  {
96  << "Patch '" << patch().name()
97  << "' corresponding to a thermal mass is not closed." << nl
98  << "Please specify the volume with the optional 'volume' entry."
99  << exit(FatalError);
100  }
101  }
102 
104 }
105 
106 
109 (
111  const fvPatch& p,
113  const fieldMapper& mapper
114 )
115 :
116  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
117  rho_(ptf.rho_),
118  Cv_(ptf.Cv_),
119  T_(ptf.T_),
120  Q_(ptf.Q_),
121  V_(ptf.V_)
122 {}
123 
124 
127 (
130 )
131 :
132  fixedValueFvPatchScalarField(ptf, iF),
133  rho_(ptf.rho_),
134  Cv_(ptf.Cv_),
135  T_(ptf.T_),
136  Q_(ptf.Q_),
137  V_(ptf.V_)
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
144 (
145  const fvPatchScalarField& ptf,
146  const fieldMapper& mapper
147 )
148 {
149  fixedValueFvPatchScalarField::map(ptf, mapper);
150 }
151 
152 
154 (
155  const fvPatchScalarField& ptf
156 )
157 {
158  fixedValueFvPatchScalarField::reset(ptf);
159 }
160 
161 
163 {
164  if
165  (
166  updated()
167  )
168  {
169  return;
170  }
171 
172  const thermophysicalTransportModel& ttm =
173  db().lookupType<thermophysicalTransportModel>
174  (
175  internalField().group()
176  );
177 
178  const scalarField Hf
179  (
180  ttm.kappaEff(patch().index())*patch().magSf()*patch().deltaCoeffs()
181  );
182 
183  const scalar Hs = rho_*Cv_*V_/db().time().deltaTValue();
184 
185  T_.value() =
186  (
187  Q_->value(db().time().value())
188  + gSum(Hf*patchInternalField())
189  + Hs*T_.oldTime().value()
190  )/(Hs + gSum(Hf));
191 
192  operator==(T_.value());
193 
194  fixedValueFvPatchScalarField::updateCoeffs();
195 }
196 
197 
199 (
200  Ostream& os
201 ) const
202 {
204  writeEntry(os, "rho", rho_);
205  writeEntry(os, "Cv", Cv_);
206  writeEntry(os, "T", T_.value());
207  writeEntry(os, Q_());
208  writeEntry(os, "volume", V_);
209  writeEntry(os, "value", *this);
210 }
211 
212 
213 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
214 
215 namespace Foam
216 {
218  (
221  );
222 }
223 
224 // ************************************************************************* //
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:91
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:237
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:257
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
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionSet dimPower
messageStream Info
const dimensionSet dimTemperature
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
static const char nl
Definition: Ostream.H:267
dictionary dict
volScalarField & p