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-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
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 
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 
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 {
113  fvPatchVectorField::operator=(Zero);
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 
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 
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 
264  // Add this side (p*area)
265  forAll(cyclicFaceCells, facei)
266  {
267  valueDiff +=p[cyclicFaceCells[facei]]*mag(initCyclicSf_[facei]);
268  }
269 
270  // Remove other side
271  forAll(nbrFaceCells, facei)
272  {
273  valueDiff -=p[nbrFaceCells[facei]]*mag(initCyclicSf_[facei]);
274  }
275 
276  if (!fBased_) //pressure based then weighted by area
277  {
278  valueDiff = valueDiff/gSum(patch().magSf());
279  }
280 
281  reduce(valueDiff, sumOp<scalar>());
282 
283  if (Pstream::master())
284  {
285  if (fBased_)
286  {
287  Info<< "Force difference = " << valueDiff << endl;
288  }
289  else
290  {
291  Info<< "Area-averaged pressure difference = "
292  << valueDiff << endl;
293  }
294  }
295 
296  if (mag(valueDiff) > mag(minThresholdValue_) || baffleActivated_)
297  {
298  openFraction_ =
299  max
300  (
301  min
302  (
303  openFraction_
304  + min
305  (
306  this->db().time().deltaT().value()/openingTime_,
307  maxOpenFractionDelta_
308  ),
309  1 - 1e-6
310  ),
311  1e-6
312  );
313 
314  baffleActivated_ = true;
315  }
316  else
317  {
318  openFraction_ = max(min(1 - 1e-6, openFraction_), 1e-6);
319  }
320 
321  if (Pstream::master())
322  {
323  Info<< "Open fraction = " << openFraction_ << endl;
324  }
325 
326  scalar areaFraction = 0.0;
327 
328  if (opening_)
329  {
330  areaFraction = openFraction_;
331  }
332  else
333  {
334  areaFraction = 1 - openFraction_;
335  }
336 
337  if (patch().boundaryMesh().mesh().moving())
338  {
339  // All areas have been recalculated already so are consistent
340  // with the new points. Note we already only do this routine
341  // once per timestep. The alternative is to use the upToDate
342  // mechanism for regIOobject + polyMesh::upToDatePoints
343  initWallSf_ = patch().Sf();
344  initCyclicSf_ = patch().boundaryMesh()[cyclicPatchLabel_].Sf();
345  nbrCyclicSf_ = refCast<const cyclicFvPatch>
346  (
347  patch().boundaryMesh()[cyclicPatchLabel_]
348  ).neighbFvPatch().Sf();
349  }
350 
351  // Update this wall patch
352  vectorField::subField Sfw = patch().patch().faceAreas();
353  vectorField newSfw((1 - areaFraction)*initWallSf_);
354  forAll(Sfw, facei)
355  {
356  Sfw[facei] = newSfw[facei];
357  }
358  const_cast<scalarField&>(patch().magSf()) = mag(patch().Sf());
359 
360  // Update owner side of cyclic
361  const_cast<vectorField&>(cyclicPatch.Sf()) = areaFraction*initCyclicSf_;
362 
363  const_cast<scalarField&>(cyclicPatch.magSf()) = mag(cyclicPatch.Sf());
364 
365  // Update neighbour side of cyclic
366  const_cast<vectorField&>(nbrPatch.Sf()) = areaFraction*nbrCyclicSf_;
367 
368  const_cast<scalarField&>(nbrPatch.magSf()) = mag(nbrPatch.Sf());
369 
370  curTimeIndex_ = this->db().time().timeIndex();
371  }
372 
373  fixedValueFvPatchVectorField::updateCoeffs();
374 }
375 
376 
378 write(Ostream& os) const
379 {
381  os.writeEntryIfDifferent<word>("p", "p", pName_);
382  os.writeEntry("cyclicPatch", cyclicPatchName_);
383  os.writeEntry("openingTime", openingTime_);
384  os.writeEntry("maxOpenFractionDelta", maxOpenFractionDelta_);
385  os.writeEntry("openFraction", openFraction_);
386  os.writeEntry("minThresholdValue", minThresholdValue_);
387  os.writeEntry("forceBased", fBased_);
388  os.writeEntry("opening", opening_);
389  writeEntry("value", os);
390 }
391 
392 
393 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
394 
395 namespace Foam
396 {
398  (
401  );
402 }
403 
404 
405 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::fvPatch::Sf
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:150
Foam::Ostream::writeEntryIfDifferent
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:248
p
volScalarField & p
Definition: createFieldRefs.H:8
cyclicFvPatch.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
surfaceFields.H
Foam::surfaceFields.
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::activePressureForceBaffleVelocityFvPatchVectorField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: activePressureForceBaffleVelocityFvPatchVectorField.C:182
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::SubField
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: Field.H:64
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::activePressureForceBaffleVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: activePressureForceBaffleVelocityFvPatchVectorField.C:378
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::polyPatch::faceCells
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:371
Foam::activePressureForceBaffleVelocityFvPatchVectorField::activePressureForceBaffleVelocityFvPatchVectorField
activePressureForceBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: activePressureForceBaffleVelocityFvPatchVectorField.C:39
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::activePressureForceBaffleVelocityFvPatchVectorField::rmap
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: activePressureForceBaffleVelocityFvPatchVectorField.C:220
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:161
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
timeIndex
label timeIndex
Definition: getTimeIndex.H:30
activePressureForceBaffleVelocityFvPatchVectorField.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::fvPatch::magSf
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:156
Foam::activePressureForceBaffleVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: activePressureForceBaffleVelocityFvPatchVectorField.C:242
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::activePressureForceBaffleVelocityFvPatchVectorField
This boundary condition is applied to the flow velocity, to simulate the opening or closure of a baff...
Definition: activePressureForceBaffleVelocityFvPatchVectorField.H:179
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54