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-2022 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."
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] =
155 (
156 "mapField",
158 );
159
160 mapTimes_[outleti] =
162 (
163 "mapTime",
165 );
166
167 ++outleti;
168 }
169
170
171 if (dict.found("value"))
172 {
174 (
175 Field<Type>("value", dict, p.size())
176 );
177 }
178 else
179 {
181 }
182}
183
184
185template<class Type>
188(
190 const fvPatch& p,
192 const fvPatchFieldMapper& mapper
193)
194:
195 fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
196 uniformValuePtr_(ptf.uniformValuePtr_.clone(p.patch())),
197 outletNames_(ptf.outletNames_),
198 offsets_(ptf.offsets_),
199 fractions_(ptf.fractions_),
200 timeDelays_(ptf.timeDelays_),
201 mapFields_(ptf.mapFields_),
202 mapTimes_(ptf.mapTimes_),
203 phiName_(ptf.phiName_),
204 curTimeIndex_(-1)
205{
206 if (mapper.direct() && !mapper.hasUnmapped())
207 {
208 // Use mapping instead of re-evaluation
209 this->map(ptf, mapper);
210 }
211 else
212 {
214 }
215}
216
217
218template<class Type>
221(
223)
224:
225 fixedValueFvPatchField<Type>(ptf),
226 uniformValuePtr_(ptf.uniformValuePtr_.clone(this->patch().patch())),
227 outletNames_(ptf.outletNames_),
228 offsets_(ptf.offsets_),
229 fractions_(ptf.fractions_),
230 timeDelays_(ptf.timeDelays_),
231 mapFields_(ptf.mapFields_),
232 mapTimes_(ptf.mapTimes_),
233 phiName_(ptf.phiName_),
234 curTimeIndex_(-1)
235{}
236
237
238template<class Type>
241(
244)
245:
246 fixedValueFvPatchField<Type>(ptf, iF),
247 uniformValuePtr_(ptf.uniformValuePtr_.clone(this->patch().patch())),
248 outletNames_(ptf.outletNames_),
249 offsets_(ptf.offsets_),
250 fractions_(ptf.fractions_),
251 timeDelays_(ptf.timeDelays_),
252 mapFields_(ptf.mapFields_),
253 mapTimes_(ptf.mapTimes_),
254 phiName_(ptf.phiName_),
255 curTimeIndex_(-1)
256{}
257
258
259// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
260
261template<class Type>
263(
264 const fvPatchFieldMapper& m
265)
266{
268
269 if (uniformValuePtr_)
270 {
271 uniformValuePtr_->autoMap(m);
272 }
273}
274
275
276template<class Type>
278(
279 const fvPatchField<Type>& ptf,
280 const labelList& addr
281)
282{
284
285 const auto& tiptf =
286 refCast<const outletMappedUniformInletFvPatchField>(ptf);
287
288 if (uniformValuePtr_)
289 {
290 uniformValuePtr_->rmap(tiptf.uniformValuePtr_(), addr);
291 }
292}
293
294
295template<class Type>
297{
298 if (this->updated())
299 {
300 return;
301 }
302
303 if (curTimeIndex_ != this->db().time().timeIndex())
304 {
305 const scalar t = this->db().time().timeOutputValue();
306
308 (
310 (
311 this->internalField()
312 )
313 );
314
315 const fvPatch& p = this->patch();
316
317 forAll(outletNames_, i)
318 {
319 const word& outletName = outletNames_[i];
320 const label outletID =
321 p.patch().boundaryMesh().findPatchID(outletName);
322
323 if (outletID < 0)
324 {
326 << "Unable to find outlet patch " << outletName
327 << abort(FatalError);
328 }
329
330
331 // Collect the map time for this outlet patch
332 DynamicList<scalar>& mapTime = mapTimes_[i];
333 scalar timeDelay = 0;
334 if (timeDelays_.set(i))
335 {
336 timeDelay = max(timeDelays_[i].value(t), scalar(0));
337 }
338 mapTime.append(t + timeDelay);
339
340
341 // Collect the map field for this outlet patch and map time
342 const fvPatchField<Type>& outletFld = f.boundaryField()[outletID];
343 DynamicList<Type>& mapField = mapFields_[i];
344
345 const auto& phi =
346 this->db().objectRegistry::template
347 lookupObject<surfaceScalarField>(phiName_);
348 const scalarField& outletPhi = phi.boundaryField()[outletID];
349 const scalar sumOutletPhi = gSum(outletPhi);
350
351 if (sumOutletPhi > SMALL)
352 {
353 Type offset(Zero);
354 if (offsets_.set(i))
355 {
356 offset = offsets_[i].value(t);
357 }
358
359 scalar fraction = 1;
360 if (fractions_.set(i))
361 {
362 fraction = fractions_[i].value(t);
363 }
364
365 mapField.append
366 (
367 gSum(outletPhi*outletFld)/sumOutletPhi*fraction
368 + offset
369 );
370 }
371 else
372 {
373 const fvPatch& outlet = p.boundaryMesh()[outletID];
374
375 mapField.append
376 (
377 gSum(outlet.magSf()*outletFld)/gSum(outlet.magSf())
378 );
379 }
380 }
381
382
383 // Map the stored fields onto inlet if the time condition is met
384 Type inletFld(Zero);
385 forAll(outletNames_, i)
386 {
387 DynamicList<scalar>& mapTime = mapTimes_[i];
388 DynamicList<Type>& mapField = mapFields_[i];
389
390 if (!mapTime.empty())
391 {
392 if (t >= mapTime.first())
393 {
394 inletFld += interpolateXY(t, mapTime, mapField);
395
396 // Remove any stored fields and times if possible
397 int i = 0;
398 while (!mapTime.empty() && t >= mapTime[i])
399 {
400 mapTime.remove(i);
401 mapField.remove(i);
402 ++i;
403 }
404 }
405 }
406 }
407
408
409 if (uniformValuePtr_)
410 {
411 this->operator==(inletFld + uniformValuePtr_->value(t));
412 }
413 else
414 {
415 this->operator==(inletFld);
416 }
417
418
419 curTimeIndex_ = this->db().time().timeIndex();
420 }
421
422
424}
425
426
427template<class Type>
429{
431
432 if (uniformValuePtr_)
433 {
434 uniformValuePtr_->writeData(os);
435 }
436 os.beginBlock("outlets");
437 forAll(outletNames_, i)
438 {
439 os.beginBlock(outletNames_[i]);
440 if (offsets_.set(i))
441 {
442 offsets_[i].writeData(os);
443 }
444 if (fractions_.set(i))
445 {
446 fractions_[i].writeData(os);
447 }
448 if (timeDelays_.set(i))
449 {
450 timeDelays_[i].writeData(os);
451 }
452 if (!mapFields_.empty())
453 {
454 mapFields_[i].writeEntry("mapField", os);
455 }
456 if (!mapTimes_.empty())
457 {
458 mapTimes_[i].writeEntry("mapTime", os);
459 }
460 os.endBlock();
461 }
462 os.endBlock();
463 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
464
465 this->writeEntry("value", os);
466}
467
468
469// ************************************************************************* //
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
T remove()
Remove and return the last element. Fatal on an empty list.
Definition: DynamicListI.H:655
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
virtual bool direct() const =0
Is it a direct (non-interpolating) mapper?
virtual bool hasUnmapped() const =0
Any unmapped values?
Generic templated field type.
Definition: Field.H:82
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:240
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:105
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:87
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:251
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
T & first()
Return the first element of the list.
Definition: UListI.H:202
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual bool write()
Write the output fields.
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...
Definition: fvPatchField.H:82
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:237
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:210
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:408
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:156
The outletMappedUniformInlet is an inlet boundary condition that.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
Interpolates y values from one curve to another with a different x distribution.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Type gSum(const FieldField< Field, Type > &f)
Field< Type > interpolateXY(const scalarField &xNew, const scalarField &xOld, const Field< Type > &yOld)
Definition: interpolateXY.C:41
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
label timeIndex
Definition: getTimeIndex.H:30
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.