activePressureForceBaffleVelocityFvPatchVectorField.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) 2016-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
31#include "volFields.H"
32#include "surfaceFields.H"
33#include "cyclicFvPatch.H"
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
37Foam::activePressureForceBaffleVelocityFvPatchVectorField::
38activePressureForceBaffleVelocityFvPatchVectorField
39(
40 const fvPatch& p,
42)
43:
44 fixedValueFvPatchVectorField(p, iF),
45 pName_("p"),
46 cyclicPatchName_(),
47 cyclicPatchLabel_(-1),
48 initWallSf_(0),
49 initCyclicSf_(0),
50 nbrCyclicSf_(0),
51 openFraction_(0),
52 openingTime_(0),
53 maxOpenFractionDelta_(0),
54 curTimeIndex_(-1),
55 minThresholdValue_(0),
56 fBased_(1),
57 baffleActivated_(0),
58 opening_(1)
59{}
60
61
62Foam::activePressureForceBaffleVelocityFvPatchVectorField::
63activePressureForceBaffleVelocityFvPatchVectorField
64(
66 const fvPatch& p,
68 const fvPatchFieldMapper& mapper
69)
70:
71 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
72 pName_(ptf.pName_),
73 cyclicPatchName_(ptf.cyclicPatchName_),
74 cyclicPatchLabel_(ptf.cyclicPatchLabel_),
75 initWallSf_(ptf.initWallSf_),
76 initCyclicSf_(ptf.initCyclicSf_),
77 nbrCyclicSf_(ptf.nbrCyclicSf_),
78 openFraction_(ptf.openFraction_),
79 openingTime_(ptf.openingTime_),
80 maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
81 curTimeIndex_(-1),
82 minThresholdValue_(ptf.minThresholdValue_),
83 fBased_(ptf.fBased_),
84 baffleActivated_(ptf.baffleActivated_),
85 opening_(ptf.opening_)
86{}
87
88
89Foam::activePressureForceBaffleVelocityFvPatchVectorField::
90activePressureForceBaffleVelocityFvPatchVectorField
91(
92 const fvPatch& p,
94 const dictionary& dict
95)
96:
97 fixedValueFvPatchVectorField(p, iF, dict, false),
98 pName_(dict.getOrDefault<word>("p", "p")),
99 cyclicPatchName_(dict.lookup("cyclicPatch")),
100 cyclicPatchLabel_(p.patch().boundaryMesh().findPatchID(cyclicPatchName_)),
101 initWallSf_(0),
102 initCyclicSf_(0),
103 nbrCyclicSf_(0),
104 openFraction_(dict.get<scalar>("openFraction")),
105 openingTime_(dict.get<scalar>("openingTime")),
106 maxOpenFractionDelta_(dict.get<scalar>("maxOpenFractionDelta")),
107 curTimeIndex_(-1),
108 minThresholdValue_(dict.get<scalar>("minThresholdValue")),
109 fBased_(dict.get<bool>("forceBased")),
110 baffleActivated_(0),
111 opening_(dict.get<bool>("opening"))
112{
114
115 if (p.size() > 0)
116 {
117 initWallSf_ = p.Sf();
118 initCyclicSf_ = p.boundaryMesh()[cyclicPatchLabel_].Sf();
119 nbrCyclicSf_ = refCast<const cyclicFvPatch>
120 (
121 p.boundaryMesh()[cyclicPatchLabel_],
122 dict
123 ).neighbFvPatch().Sf();
124 }
125
126 dict.readIfPresent("p", pName_);
127}
128
129
130Foam::activePressureForceBaffleVelocityFvPatchVectorField::
131activePressureForceBaffleVelocityFvPatchVectorField
132(
134)
135:
136 fixedValueFvPatchVectorField(ptf),
137 pName_(ptf.pName_),
138 cyclicPatchName_(ptf.cyclicPatchName_),
139 cyclicPatchLabel_(ptf.cyclicPatchLabel_),
140 initWallSf_(ptf.initWallSf_),
141 initCyclicSf_(ptf.initCyclicSf_),
142 nbrCyclicSf_(ptf.nbrCyclicSf_),
143 openFraction_(ptf.openFraction_),
144 openingTime_(ptf.openingTime_),
145 maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
146 curTimeIndex_(-1),
147 minThresholdValue_(ptf.minThresholdValue_),
148 fBased_(ptf.fBased_),
149 baffleActivated_(ptf.baffleActivated_),
150 opening_(ptf.opening_)
151{}
152
153
154Foam::activePressureForceBaffleVelocityFvPatchVectorField::
155activePressureForceBaffleVelocityFvPatchVectorField
156(
159)
160:
161 fixedValueFvPatchVectorField(ptf, iF),
162 pName_(ptf.pName_),
163 cyclicPatchName_(ptf.cyclicPatchName_),
164 cyclicPatchLabel_(ptf.cyclicPatchLabel_),
165 initWallSf_(ptf.initWallSf_),
166 initCyclicSf_(ptf.initCyclicSf_),
167 nbrCyclicSf_(ptf.nbrCyclicSf_),
168 openFraction_(ptf.openFraction_),
169 openingTime_(ptf.openingTime_),
170 maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
171 curTimeIndex_(-1),
172 minThresholdValue_(ptf.minThresholdValue_),
173 fBased_(ptf.fBased_),
174 baffleActivated_(ptf.baffleActivated_),
175 opening_(ptf.opening_)
176{}
177
178
179// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
180
182(
183 const fvPatchFieldMapper& m
184)
185{
186 fixedValueFvPatchVectorField::autoMap(m);
187
188 //- Note: cannot map field from cyclic patch anyway so just recalculate
189 // Areas should be consistent when doing autoMap except in case of
190 // topo changes.
191 //- Note: we don't want to use Sf here since triggers rebuilding of
192 // fvMesh::S() which will give problems when mapped (since already
193 // on new mesh)
194 forAll(patch().boundaryMesh().mesh().faceAreas(), i)
195 {
196 if (mag(patch().boundaryMesh().mesh().faceAreas()[i]) == 0)
197 {
198 Info << "faceArea[active] "<< i << endl;
199 }
200 }
201
202 if (patch().size() > 0)
203 {
204 const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
205 initWallSf_ = patch().patchSlice(areas);
206 initCyclicSf_ =
207 patch().boundaryMesh()[cyclicPatchLabel_].patchSlice(areas);
208 nbrCyclicSf_ = refCast<const cyclicFvPatch>
209 (
210 patch().boundaryMesh()
211 [
212 cyclicPatchLabel_
213 ]
214 ).neighbFvPatch().patch().patchSlice(areas);
215 }
216}
217
218
220(
221 const fvPatchVectorField& ptf,
222 const labelList& addr
223)
224{
225 fixedValueFvPatchVectorField::rmap(ptf, addr);
226
227 // See autoMap.
228 const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
229 initWallSf_ = patch().patchSlice(areas);
230 initCyclicSf_ =
231 patch().boundaryMesh()[cyclicPatchLabel_].patchSlice(areas);
232 nbrCyclicSf_ = refCast<const cyclicFvPatch>
233 (
234 patch().boundaryMesh()
235 [
236 cyclicPatchLabel_
237 ]
238 ).neighbFvPatch().patch().patchSlice(areas);
239}
240
241
243{
244 if (updated())
245 {
246 return;
247 }
248
249 // Execute the change to the openFraction only once per time-step
250 if (curTimeIndex_ != this->db().time().timeIndex())
251 {
252 const volScalarField& p =
253 db().lookupObject<volScalarField>(pName_);
254
255 const fvPatch& cyclicPatch = patch().boundaryMesh()[cyclicPatchLabel_];
256 const labelList& cyclicFaceCells = cyclicPatch.patch().faceCells();
257 const fvPatch& nbrPatch =
258 refCast<const cyclicFvPatch>(cyclicPatch).neighbFvPatch();
259
260 const labelList& nbrFaceCells = nbrPatch.patch().faceCells();
261
262 scalar valueDiff = 0;
263 scalar area = 0;
264
265 // Add this side (p*area)
266 forAll(cyclicFaceCells, facei)
267 {
268 valueDiff +=p[cyclicFaceCells[facei]]*mag(initCyclicSf_[facei]);
269 area += mag(initCyclicSf_[facei]);
270 }
271
272 // Remove other side
273 forAll(nbrFaceCells, facei)
274 {
275 valueDiff -=p[nbrFaceCells[facei]]*mag(initCyclicSf_[facei]);
276 }
277
278 if (!fBased_) //pressure based then weighted by area
279 {
280 valueDiff = valueDiff/(area + VSMALL);
281 }
282
283 reduce(valueDiff, sumOp<scalar>());
284
285 if (Pstream::master())
286 {
287 if (fBased_)
288 {
289 Info<< "Force difference (threshold) = " << valueDiff
290 << "(" << minThresholdValue_ << ")" << endl;
291 }
292 else
293 {
294 Info<< "Area-averaged pressure difference (threshold) = "
295 << valueDiff << "(" << minThresholdValue_ << ")" << endl;
296 }
297 }
298
299 if (mag(valueDiff) > mag(minThresholdValue_) || baffleActivated_)
300 {
301 openFraction_ =
302 max
303 (
304 min
305 (
306 openFraction_
307 + min
308 (
309 this->db().time().deltaT().value()/openingTime_,
310 maxOpenFractionDelta_
311 ),
312 1 - 1e-6
313 ),
314 1e-6
315 );
316
317 baffleActivated_ = true;
318 }
319 else
320 {
321 openFraction_ = max(min(1 - 1e-6, openFraction_), 1e-6);
322 }
323
324 if (Pstream::master())
325 {
326 Info<< "Open fraction = " << openFraction_ << endl;
327 }
328
329 scalar areaFraction = 0.0;
330
331 if (opening_)
332 {
333 areaFraction = openFraction_;
334 }
335 else
336 {
337 areaFraction = 1 - openFraction_;
338 }
339
340 if (patch().boundaryMesh().mesh().moving())
341 {
342 // All areas have been recalculated already so are consistent
343 // with the new points. Note we already only do this routine
344 // once per timestep. The alternative is to use the upToDate
345 // mechanism for regIOobject + polyMesh::upToDatePoints
346 initWallSf_ = patch().Sf();
347 initCyclicSf_ = patch().boundaryMesh()[cyclicPatchLabel_].Sf();
348 nbrCyclicSf_ = refCast<const cyclicFvPatch>
349 (
350 patch().boundaryMesh()[cyclicPatchLabel_]
351 ).neighbFvPatch().Sf();
352 }
353
354 // Update this wall patch
355 vectorField::subField Sfw = patch().patch().faceAreas();
356 vectorField newSfw((1 - areaFraction)*initWallSf_);
357 forAll(Sfw, facei)
358 {
359 Sfw[facei] = newSfw[facei];
360 }
361 const_cast<scalarField&>(patch().magSf()) = mag(patch().Sf());
362
363 // Update owner side of cyclic
364 const_cast<vectorField&>(cyclicPatch.Sf()) = areaFraction*initCyclicSf_;
365
366 const_cast<scalarField&>(cyclicPatch.magSf()) = mag(cyclicPatch.Sf());
367
368 // Update neighbour side of cyclic
369 const_cast<vectorField&>(nbrPatch.Sf()) = areaFraction*nbrCyclicSf_;
370
371 const_cast<scalarField&>(nbrPatch.magSf()) = mag(nbrPatch.Sf());
372
373 curTimeIndex_ = this->db().time().timeIndex();
374 }
375
376 fixedValueFvPatchVectorField::updateCoeffs();
377}
378
379
381write(Ostream& os) const
382{
384 os.writeEntryIfDifferent<word>("p", "p", pName_);
385 os.writeEntry("cyclicPatch", cyclicPatchName_);
386 os.writeEntry("openingTime", openingTime_);
387 os.writeEntry("maxOpenFractionDelta", maxOpenFractionDelta_);
388 os.writeEntry("openFraction", openFraction_);
389 os.writeEntry("minThresholdValue", minThresholdValue_);
390 os.writeEntry("forceBased", fBased_);
391 os.writeEntry("opening", opening_);
392 writeEntry("value", os);
393}
394
395
396// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
397
398namespace Foam
399{
401 (
404 );
405}
406
407
408// ************************************************************************* //
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: SubField.H:62
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
This boundary condition is applied to the flow velocity, to simulate the opening or closure of a baff...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< vector > &)
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
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:167
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:150
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:371
Lookup type of boundary radiation properties.
Definition: lookup.H:66
splitCell * master() const
Definition: splitCell.H:113
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
bool
Definition: EEqn.H:20
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.