Loading...
Searching...
No Matches
outletMappedUniformInletFvPatchField.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "volFields.H"
31#include "surfaceFields.H"
32#include "interpolateXY.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class Type>
39(
40 const fvPatch& p,
42)
43:
44 fixedValueFvPatchField<Type>(p, iF),
45 uniformValuePtr_(nullptr),
46 outletNames_(),
47 offsets_(),
48 fractions_(),
49 timeDelays_(),
50 mapFields_(),
51 mapTimes_(),
52 phiName_("phi"),
53 curTimeIndex_(-1)
54{}
55
56
57template<class Type>
60(
61 const fvPatch& p,
63 const dictionary& dict
64)
65:
66 fixedValueFvPatchField<Type>(p, iF, dict),
67 uniformValuePtr_
68 (
69 PatchFunction1<Type>::NewIfPresent
70 (
71 p.patch(),
72 "uniformValue",
73 dict
74 )
75 ),
76 outletNames_(),
77 offsets_(),
78 fractions_(),
79 timeDelays_(),
80 mapFields_(),
81 mapTimes_(),
82 phiName_(dict.getOrDefault<word>("phi", "phi")),
83 curTimeIndex_(-1)
84{
85 const dictionary& outletDict = dict.subDict("outlets");
86
87 if (outletDict.empty())
88 {
89 FatalIOErrorInFunction(outletDict)
90 << "outlets dictionary is empty."
91 << exit(FatalIOError);
92 }
93
94 outletNames_.setSize(outletDict.size());
95 offsets_.setSize(outletDict.size());
96 fractions_.setSize(outletDict.size());
97 timeDelays_.setSize(outletDict.size());
98 mapFields_.setSize(outletDict.size());
99 mapTimes_.setSize(outletDict.size());
100
101 label outleti = 0;
102 for (const entry& dEntry : outletDict)
103 {
104 const word& key = dEntry.keyword();
105
106 if (!dEntry.isDict())
107 {
108 FatalIOErrorInFunction(outletDict)
109 << "Entry " << key << " is not a dictionary." << nl
110 << exit(FatalIOError);
111 }
112
113 const dictionary& subDict = dEntry.dict();
114
115 outletNames_[outleti] = key;
116
117 offsets_.set
118 (
119 outleti,
121 (
122 "offset",
123 subDict,
125 &this->db()
126 )
127 );
128
129 fractions_.set
130 (
131 outleti,
133 (
134 "fraction",
135 subDict,
137 &this->db()
138 )
139 );
140
141 timeDelays_.set
142 (
143 outleti,
145 (
146 "timeDelay",
147 subDict,
149 &this->db()
150 )
151 );
152
153 mapFields_[outleti] =
154 subDict.getOrDefault<DynamicList<Type>>
155 (
156 "mapField",
157 DynamicList<Type>()
158 );
159
160 mapTimes_[outleti] =
161 subDict.getOrDefault<DynamicList<scalar>>
162 (
163 "mapTime",
164 DynamicList<scalar>()
165 );
166
167 ++outleti;
168 }
169
170
171 if (!this->readValueEntry(dict))
172 {
173 // Fallback: set to the internal field
174 this->extrapolateInternal();
175 }
176}
177
178
179template<class Type>
182(
183 const outletMappedUniformInletFvPatchField<Type>& ptf,
184 const fvPatch& p,
185 const DimensionedField<Type, volMesh>& iF,
186 const fvPatchFieldMapper& mapper
187)
188:
189 fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
190 uniformValuePtr_(ptf.uniformValuePtr_.clone(p.patch())),
191 outletNames_(ptf.outletNames_),
192 offsets_(ptf.offsets_),
193 fractions_(ptf.fractions_),
194 timeDelays_(ptf.timeDelays_),
195 mapFields_(ptf.mapFields_),
196 mapTimes_(ptf.mapTimes_),
197 phiName_(ptf.phiName_),
198 curTimeIndex_(-1)
199{
200 if (mapper.direct() && !mapper.hasUnmapped())
201 {
202 // Use mapping instead of re-evaluation
203 this->map(ptf, mapper);
204 }
205 else
206 {
207 // Fallback: set to the internal field
208 this->extrapolateInternal();
209 }
210}
211
212
213template<class Type>
216(
218)
219:
220 fixedValueFvPatchField<Type>(ptf),
221 uniformValuePtr_(ptf.uniformValuePtr_.clone(this->patch().patch())),
222 outletNames_(ptf.outletNames_),
223 offsets_(ptf.offsets_),
224 fractions_(ptf.fractions_),
225 timeDelays_(ptf.timeDelays_),
226 mapFields_(ptf.mapFields_),
227 mapTimes_(ptf.mapTimes_),
228 phiName_(ptf.phiName_),
229 curTimeIndex_(-1)
230{}
231
232
233template<class Type>
236(
239)
240:
241 fixedValueFvPatchField<Type>(ptf, iF),
242 uniformValuePtr_(ptf.uniformValuePtr_.clone(this->patch().patch())),
243 outletNames_(ptf.outletNames_),
244 offsets_(ptf.offsets_),
245 fractions_(ptf.fractions_),
246 timeDelays_(ptf.timeDelays_),
247 mapFields_(ptf.mapFields_),
248 mapTimes_(ptf.mapTimes_),
249 phiName_(ptf.phiName_),
250 curTimeIndex_(-1)
251{}
252
253
254// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
255
256template<class Type>
258(
259 const fvPatchFieldMapper& m
260)
261{
263
264 if (uniformValuePtr_)
266 uniformValuePtr_->autoMap(m);
267 }
268}
269
270
271template<class Type>
273(
274 const fvPatchField<Type>& ptf,
275 const labelList& addr
276)
277{
279
280 const auto& tiptf =
282
283 if (uniformValuePtr_)
285 uniformValuePtr_->rmap(tiptf.uniformValuePtr_(), addr);
286 }
287}
288
289
290template<class Type>
292{
293 if (this->updated())
294 {
295 return;
296 }
297
298 if (curTimeIndex_ != this->db().time().timeIndex())
299 {
300 const scalar t = this->db().time().timeOutputValue();
301
303 (
305 (
306 this->internalField()
307 )
308 );
309
310 const fvPatch& p = this->patch();
311
312 forAll(outletNames_, i)
313 {
314 const word& outletName = outletNames_[i];
315 const label outletID =
316 p.patch().boundaryMesh().findPatchID(outletName);
317
318 if (outletID < 0)
319 {
321 << "Unable to find outlet patch " << outletName
322 << abort(FatalError);
323 }
324
325
326 // Collect the map time for this outlet patch
327 DynamicList<scalar>& mapTime = mapTimes_[i];
328 scalar timeDelay = 0;
329 if (timeDelays_.set(i))
330 {
331 timeDelay = max(timeDelays_[i].value(t), scalar(0));
332 }
333 mapTime.append(t + timeDelay);
334
335
336 // Collect the map field for this outlet patch and map time
337 const fvPatchField<Type>& outletFld = f.boundaryField()[outletID];
338 DynamicList<Type>& mapField = mapFields_[i];
339
340 const auto& phi =
341 this->db().objectRegistry::template
342 lookupObject<surfaceScalarField>(phiName_);
343 const scalarField& outletPhi = phi.boundaryField()[outletID];
344 const scalar sumOutletPhi = gSum(outletPhi);
345
346 if (sumOutletPhi > SMALL)
347 {
348 Type offset(Zero);
349 if (offsets_.set(i))
350 {
351 offset = offsets_[i].value(t);
352 }
353
354 scalar fraction = 1;
355 if (fractions_.set(i))
356 {
357 fraction = fractions_[i].value(t);
358 }
359
360 mapField.append
361 (
362 gSum(outletPhi*outletFld)/sumOutletPhi*fraction
363 + offset
364 );
365 }
366 else
367 {
368 const fvPatch& outlet = p.boundaryMesh()[outletID];
369
370 mapField.append
371 (
372 gWeightedAverage(outlet.magSf(), outletFld)
373 );
374 }
375 }
376
377
378 // Map the stored fields onto inlet if the time condition is met
379 Type inletFld(Zero);
380 forAll(outletNames_, i)
381 {
382 DynamicList<scalar>& mapTime = mapTimes_[i];
383 DynamicList<Type>& mapField = mapFields_[i];
384
385 if (!mapTime.empty())
386 {
387 if (t >= mapTime.first())
388 {
389 inletFld += interpolateXY(t, mapTime, mapField);
390
391 // Remove any stored fields and times if possible
392 int i = 0;
393 while (!mapTime.empty() && t >= mapTime[i])
394 {
395 mapTime.remove(i);
396 mapField.remove(i);
397 ++i;
398 }
399 }
400 }
401 }
402
403
404 if (uniformValuePtr_)
405 {
406 this->operator==(inletFld + uniformValuePtr_->value(t));
407 }
408 else
409 {
410 this->operator==(inletFld);
411 }
412
413
414 curTimeIndex_ = this->db().time().timeIndex();
415 }
417
419}
420
421
422template<class Type>
424{
426
427 if (uniformValuePtr_)
428 {
429 uniformValuePtr_->writeData(os);
430 }
431 os.beginBlock("outlets");
432 forAll(outletNames_, i)
433 {
434 os.beginBlock(outletNames_[i]);
435 if (offsets_.set(i))
436 {
437 offsets_[i].writeData(os);
438 }
439 if (fractions_.set(i))
440 {
441 fractions_[i].writeData(os);
442 }
443 if (timeDelays_.set(i))
444 {
445 timeDelays_[i].writeData(os);
446 }
447 if (!mapFields_.empty())
448 {
449 mapFields_[i].writeEntry("mapField", os);
450 }
451 if (!mapTimes_.empty())
452 {
453 mapTimes_[i].writeEntry("mapTime", os);
454 }
455 os.endBlock();
456 }
457 os.endBlock();
458 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
459
461}
462
463
464// ************************************************************************* //
for(const label curEdgei :curPointEdges)
bool empty() const noexcept
True if the list is empty.
Definition DLListBase.H:189
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void append(const T &val)
Append an element at the end of the list.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
T remove()
Remove and return the last element. Fatal on an empty list.
void append(const T &val)
Copy append an element to the end of this list.
virtual bool direct() const =0
Is it a direct (non-interpolating) mapper?
virtual bool hasUnmapped() const =0
Any unmapped values?
static autoPtr< Function1< Type > > NewIfPresent(const word &entryName, const dictionary &dict, const word &redirectType, const objectRegistry *obrPtr=nullptr)
An optional selector, with fallback redirection.
Generic GeometricField class.
void setSize(label n)
Alias for resize().
Definition List.H:536
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
virtual Ostream & endBlock()
Write end block group.
Definition Ostream.C:108
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition Ostream.H:346
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
T & first()
Access first element of the list, position [0].
Definition UList.H:957
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
fixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
void extrapolateInternal()
Assign the patch field from the internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
const scalarField & magSf() const
Return face area magnitudes, like the fvMesh::magSf() method.
Definition fvPatch.C:131
The outletMappedUniformInlet is an inlet boundary condition that.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
outletMappedUniformInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
Interpolates y values from one curve to another with a different x distribution.
const std::string patch
OpenFOAM patch number as a std::string.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
Type gSum(const FieldField< Field, Type > &f)
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
List< label > labelList
A List of labels.
Definition List.H:62
Field< Type > interpolateXY(const scalarField &xNew, const scalarField &xOld, const Field< Type > &yOld)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
errorManip< error > abort(error &err)
Definition errorManip.H:139
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label timeIndex
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.