activePressureForceBaffleVelocityFvPatchVectorField.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) 2011-2022 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 
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "cyclicFvPatch.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  fixedValueFvPatchVectorField(p, iF),
42  pName_("p"),
43  cyclicPatchName_(),
44  cyclicPatchLabel_(-1),
45  orientation_(1),
46  initWallSf_(0),
47  initCyclicSf_(0),
48  nbrCyclicSf_(0),
49  openFraction_(0),
50  openingTime_(0),
51  maxOpenFractionDelta_(0),
52  curTimeIndex_(-1),
53  minThresholdValue_(0),
54  fBased_(1),
55  baffleActivated_(0)
56 {}
57 
58 
61 (
62  const fvPatch& p,
64  const dictionary& dict
65 )
66 :
67  fixedValueFvPatchVectorField(p, iF, dict, false),
68  pName_(dict.lookupOrDefault<word>("p", "p")),
69  cyclicPatchName_(dict.lookup("cyclicPatch")),
70  cyclicPatchLabel_(p.patch().boundaryMesh().findPatchID(cyclicPatchName_)),
71  orientation_(dict.lookup<label>("orientation")),
72  initWallSf_(0),
73  initCyclicSf_(0),
74  nbrCyclicSf_(0),
75  openFraction_(dict.lookup<scalar>("openFraction")),
76  openingTime_(dict.lookup<scalar>("openingTime")),
77  maxOpenFractionDelta_(dict.lookup<scalar>("maxOpenFractionDelta")),
78  curTimeIndex_(-1),
79  minThresholdValue_(dict.lookup<scalar>("minThresholdValue")),
80  fBased_(readBool(dict.lookup("forceBased"))),
81  baffleActivated_(0)
82 {
83  fvPatchVectorField::operator=(Zero);
84 
85  if (p.size() > 0)
86  {
87  initWallSf_ = p.Sf();
88  initCyclicSf_ = p.boundaryMesh()[cyclicPatchLabel_].Sf();
89  nbrCyclicSf_ = refCast<const cyclicFvPatch>
90  (
91  p.boundaryMesh()[cyclicPatchLabel_]
92  ).neighbFvPatch().Sf();
93  }
94 
95  if (dict.found("p"))
96  {
97  dict.lookup("p") >> pName_;
98  }
99 }
100 
101 
104 (
106  const fvPatch& p,
108  const fvPatchFieldMapper& mapper
109 )
110 :
111  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
112  pName_(ptf.pName_),
113  cyclicPatchName_(ptf.cyclicPatchName_),
114  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
115  orientation_(ptf.orientation_),
116  initWallSf_(ptf.initWallSf_),
117  initCyclicSf_(ptf.initCyclicSf_),
118  nbrCyclicSf_(ptf.nbrCyclicSf_),
119  openFraction_(ptf.openFraction_),
120  openingTime_(ptf.openingTime_),
121  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
122  curTimeIndex_(-1),
123  minThresholdValue_(ptf.minThresholdValue_),
124  fBased_(ptf.fBased_),
125  baffleActivated_(ptf.baffleActivated_)
126 {}
127 
128 
131 (
134 )
135 :
136  fixedValueFvPatchVectorField(ptf, iF),
137  pName_(ptf.pName_),
138  cyclicPatchName_(ptf.cyclicPatchName_),
139  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
140  orientation_(ptf.orientation_),
141  initWallSf_(ptf.initWallSf_),
142  initCyclicSf_(ptf.initCyclicSf_),
143  nbrCyclicSf_(ptf.nbrCyclicSf_),
144  openFraction_(ptf.openFraction_),
145  openingTime_(ptf.openingTime_),
146  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
147  curTimeIndex_(-1),
148  minThresholdValue_(ptf.minThresholdValue_),
149  fBased_(ptf.fBased_),
150  baffleActivated_(ptf.baffleActivated_)
151 {}
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
157 (
158  const fvPatchFieldMapper& m
159 )
160 {
161  fixedValueFvPatchVectorField::autoMap(m);
162 
163  //- Note: cannot map field from cyclic patch anyway so just recalculate
164  // Areas should be consistent when doing autoMap except in case of
165  // topo changes.
166  //- Note: we don't want to use Sf here since triggers rebuilding of
167  // fvMesh::S() which will give problems when mapped (since already
168  // on new mesh)
169  forAll(patch().boundaryMesh().mesh().faceAreas(), i)
170  {
171  if (mag(patch().boundaryMesh().mesh().faceAreas()[i]) == 0)
172  {
173  Info << "faceArea[active] "<< i << endl;
174  }
175  }
176  if (patch().size() > 0)
177  {
178  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
179  initWallSf_ = patch().patchSlice(areas);
180  initCyclicSf_ = patch().boundaryMesh()
181  [
182  cyclicPatchLabel_
183  ].patchSlice(areas);
184  nbrCyclicSf_ = refCast<const cyclicFvPatch>
185  (
186  patch().boundaryMesh()
187  [
188  cyclicPatchLabel_
189  ]
190  ).neighbFvPatch().patch().patchSlice(areas);
191  }
192 }
193 
194 
196 (
197  const fvPatchVectorField& ptf,
198  const labelList& addr
199 )
200 {
201  fixedValueFvPatchVectorField::rmap(ptf, addr);
202 
203  // See autoMap.
204  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
205  initWallSf_ = patch().patchSlice(areas);
206  initCyclicSf_ = patch().boundaryMesh()
207  [
208  cyclicPatchLabel_
209  ].patchSlice(areas);
210  nbrCyclicSf_ = refCast<const cyclicFvPatch>
211  (
212  patch().boundaryMesh()
213  [
214  cyclicPatchLabel_
215  ]
216  ).neighbFvPatch().patch().patchSlice(areas);
217 }
218 
219 
221 (
222  const fvPatchVectorField& ptf
223 )
224 {
225  fixedValueFvPatchVectorField::reset(ptf);
226 
227  // See autoMap.
228  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
229  initWallSf_ = patch().patchSlice(areas);
230  initCyclicSf_ = patch().boundaryMesh()
231  [
232  cyclicPatchLabel_
233  ].patchSlice(areas);
234  nbrCyclicSf_ = refCast<const cyclicFvPatch>
235  (
236  patch().boundaryMesh()
237  [
238  cyclicPatchLabel_
239  ]
240  ).neighbFvPatch().patch().patchSlice(areas);
241 }
242 
243 
245 {
246  if (updated())
247  {
248  return;
249  }
250  // Execute the change to the openFraction only once per time-step
251  if (curTimeIndex_ != this->db().time().timeIndex())
252  {
253  const volScalarField& p = db().lookupObject<volScalarField>
254  (
255  pName_
256  );
257 
258  const fvPatch& cyclicPatch = patch().boundaryMesh()[cyclicPatchLabel_];
259  const labelList& cyclicFaceCells = cyclicPatch.patch().faceCells();
260  const fvPatch& nbrPatch = refCast<const cyclicFvPatch>
261  (
262  cyclicPatch
263  ).neighbFvPatch();
264 
265  const labelList& nbrFaceCells = nbrPatch.patch().faceCells();
266 
267  scalar valueDiff = 0;
268 
269  if (fBased_)
270  {
271  // Add this side
272  forAll(cyclicFaceCells, facei)
273  {
274  valueDiff +=p[cyclicFaceCells[facei]]*mag(initCyclicSf_[facei]);
275  }
276 
277  // Remove other side
278  forAll(nbrFaceCells, facei)
279  {
280  valueDiff -=p[nbrFaceCells[facei]]*mag(initCyclicSf_[facei]);
281  }
282 
283  Info<< "Force difference = " << valueDiff << endl;
284  }
285  else // pressure based
286  {
287  forAll(cyclicFaceCells, facei)
288  {
289  valueDiff += p[cyclicFaceCells[facei]];
290  }
291 
292  forAll(nbrFaceCells, facei)
293  {
294  valueDiff -= p[nbrFaceCells[facei]];
295  }
296 
297  Info<< "Pressure difference = " << valueDiff << endl;
298  }
299 
300  if ((mag(valueDiff) > mag(minThresholdValue_)) || baffleActivated_)
301  {
302  openFraction_ =
303  max(
304  min(
305  openFraction_
306  + min
307  (
308  this->db().time().deltaT().value()/openingTime_,
309  maxOpenFractionDelta_
310  )*(orientation_),
311  1 - 1e-6
312  ),
313  1e-6
314  );
315 
316  baffleActivated_ = true;
317  }
318  else
319  {
320  openFraction_ = max(min(1 - 1e-6, openFraction_), 1e-6);
321  }
322 
323  Info<< "Open fraction = " << openFraction_ << endl;
324 
325  vectorField::subField Sfw = patch().patch().faceAreas();
326  vectorField newSfw((1 - openFraction_)*initWallSf_);
327  forAll(Sfw, facei)
328  {
329  Sfw[facei] = newSfw[facei];
330  }
331  const_cast<scalarField&>(patch().magSf()) = mag(patch().Sf());
332 
333  // Update owner side of cyclic
334  const_cast<vectorField&>(cyclicPatch.Sf()) =
335  openFraction_*initCyclicSf_;
336 
337  const_cast<scalarField&>(cyclicPatch.magSf()) =
338  mag(cyclicPatch.Sf());
339 
340  // Update neighbour side of cyclic
341  const_cast<vectorField&>(nbrPatch.Sf()) =
342  openFraction_*nbrCyclicSf_;
343 
344  const_cast<scalarField&>(nbrPatch.magSf()) =
345  mag(nbrPatch.Sf());
346 
347  curTimeIndex_ = this->db().time().timeIndex();
348  }
349 
350  fixedValueFvPatchVectorField::updateCoeffs();
351 }
352 
353 
355 write(Ostream& os) const
356 {
358  writeEntryIfDifferent<word>(os, "p", "p", pName_);
359  writeEntry(os, "cyclicPatch", cyclicPatchName_);
360  writeEntry(os, "orientation", orientation_);
361  writeEntry(os, "openingTime", openingTime_);
362  writeEntry(os, "maxOpenFractionDelta", maxOpenFractionDelta_);
363  writeEntry(os, "openFraction", openFraction_);
364  writeEntry(os, "minThresholdValue", minThresholdValue_);
365  writeEntry(os, "forceBased", fBased_);
366  writeEntry(os, "value", *this);
367 }
368 
369 
370 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
371 
372 namespace Foam
373 {
375  (
378  );
379 }
380 
381 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:181
This boundary condition is applied to the flow velocity, to simulate the opening or closure of a baff...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:142
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:297
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
label findPatchID(const word &patchName) const
Find patch index given a name.
bool readBool(Istream &)
Definition: boolIO.C:60
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Pre-declare related SubField type.
Definition: Field.H:60
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
A class for handling words, derived from string.
Definition: word.H:59
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:139
activePressureForceBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Foam::fvPatchFieldMapper.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
static const zero Zero
Definition: zero.H:97
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:136
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
label timeIndex
Definition: getTimeIndex.H:4
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Namespace for OpenFOAM.
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864