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-2024 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  const dictionary& dict
40 )
41 :
42  fixedValueFvPatchVectorField(p, iF, dict, false),
43  pName_(dict.lookupOrDefault<word>("p", "p")),
44  cyclicPatchName_(dict.lookup("cyclicPatch")),
45  cyclicPatchLabel_(p.patch().boundaryMesh().findIndex(cyclicPatchName_)),
46  orientation_(dict.lookup<label>("orientation")),
47  initWallSf_(0),
48  initCyclicSf_(0),
49  nbrCyclicSf_(0),
50  openFraction_(dict.lookup<scalar>("openFraction", unitFraction)),
51  openingTime_(dict.lookup<scalar>("openingTime", dimTime)),
52  maxOpenFractionDelta_
53  (
54  dict.lookup<scalar>("maxOpenFractionDelta", unitFraction)
55  ),
56  curTimeIndex_(-1),
57  minThresholdValue_(dict.lookup<scalar>("minThresholdValue", dimPressure)),
58  fBased_(dict.lookup<bool>("forceBased")),
59  baffleActivated_(0)
60 {
62 
63  if (p.size() > 0)
64  {
65  initWallSf_ = p.Sf();
66  initCyclicSf_ = p.boundaryMesh()[cyclicPatchLabel_].Sf();
67  nbrCyclicSf_ = refCast<const cyclicFvPatch>
68  (
69  p.boundaryMesh()[cyclicPatchLabel_]
70  ).neighbFvPatch().Sf();
71  }
72 
73  if (dict.found("p"))
74  {
75  dict.lookup("p") >> pName_;
76  }
77 }
78 
79 
82 (
84  const fvPatch& p,
86  const fieldMapper& mapper
87 )
88 :
89  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
90  pName_(ptf.pName_),
91  cyclicPatchName_(ptf.cyclicPatchName_),
92  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
93  orientation_(ptf.orientation_),
94  initWallSf_(ptf.initWallSf_),
95  initCyclicSf_(ptf.initCyclicSf_),
96  nbrCyclicSf_(ptf.nbrCyclicSf_),
97  openFraction_(ptf.openFraction_),
98  openingTime_(ptf.openingTime_),
99  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
100  curTimeIndex_(-1),
101  minThresholdValue_(ptf.minThresholdValue_),
102  fBased_(ptf.fBased_),
103  baffleActivated_(ptf.baffleActivated_)
104 {}
105 
106 
109 (
112 )
113 :
114  fixedValueFvPatchVectorField(ptf, iF),
115  pName_(ptf.pName_),
116  cyclicPatchName_(ptf.cyclicPatchName_),
117  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
118  orientation_(ptf.orientation_),
119  initWallSf_(ptf.initWallSf_),
120  initCyclicSf_(ptf.initCyclicSf_),
121  nbrCyclicSf_(ptf.nbrCyclicSf_),
122  openFraction_(ptf.openFraction_),
123  openingTime_(ptf.openingTime_),
124  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
125  curTimeIndex_(-1),
126  minThresholdValue_(ptf.minThresholdValue_),
127  fBased_(ptf.fBased_),
128  baffleActivated_(ptf.baffleActivated_)
129 {}
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
135 (
136  const fvPatchVectorField& ptf,
137  const fieldMapper& mapper
138 )
139 {
140  fixedValueFvPatchVectorField::map(ptf, mapper);
141 
142  //- Note: cannot map field from cyclic patch anyway so just recalculate
143  // Areas should be consistent when doing map except in case of topo
144  // changes.
145  //- Note: we don't want to use Sf here since triggers rebuilding of
146  // fvMesh::S() which will give problems when mapped (since already
147  // on new mesh)
148 
149  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
150  initWallSf_ = patch().patchSlice(areas);
151  initCyclicSf_ = patch().boundaryMesh()
152  [
153  cyclicPatchLabel_
154  ].patchSlice(areas);
155  nbrCyclicSf_ = refCast<const cyclicFvPatch>
156  (
157  patch().boundaryMesh()
158  [
159  cyclicPatchLabel_
160  ]
161  ).neighbFvPatch().patch().patchSlice(areas);
162 }
163 
164 
166 (
167  const fvPatchVectorField& ptf
168 )
169 {
170  fixedValueFvPatchVectorField::reset(ptf);
171 
172  // See rmap.
173  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
174  initWallSf_ = patch().patchSlice(areas);
175  initCyclicSf_ = patch().boundaryMesh()
176  [
177  cyclicPatchLabel_
178  ].patchSlice(areas);
179  nbrCyclicSf_ = refCast<const cyclicFvPatch>
180  (
181  patch().boundaryMesh()
182  [
183  cyclicPatchLabel_
184  ]
185  ).neighbFvPatch().patch().patchSlice(areas);
186 }
187 
188 
190 {
191  if (updated())
192  {
193  return;
194  }
195  // Execute the change to the openFraction only once per time-step
196  if (curTimeIndex_ != this->db().time().timeIndex())
197  {
198  const volScalarField& p = db().lookupObject<volScalarField>
199  (
200  pName_
201  );
202 
203  const fvPatch& cyclicPatch = patch().boundaryMesh()[cyclicPatchLabel_];
204  const labelList& cyclicFaceCells = cyclicPatch.patch().faceCells();
205  const fvPatch& nbrPatch = refCast<const cyclicFvPatch>
206  (
207  cyclicPatch
208  ).neighbFvPatch();
209 
210  const labelList& nbrFaceCells = nbrPatch.patch().faceCells();
211 
212  scalar valueDiff = 0;
213 
214  if (fBased_)
215  {
216  // Add this side
217  forAll(cyclicFaceCells, facei)
218  {
219  valueDiff +=p[cyclicFaceCells[facei]]*mag(initCyclicSf_[facei]);
220  }
221 
222  // Remove other side
223  forAll(nbrFaceCells, facei)
224  {
225  valueDiff -=p[nbrFaceCells[facei]]*mag(initCyclicSf_[facei]);
226  }
227 
228  Info<< "Force difference = " << valueDiff << endl;
229  }
230  else // pressure based
231  {
232  forAll(cyclicFaceCells, facei)
233  {
234  valueDiff += p[cyclicFaceCells[facei]];
235  }
236 
237  forAll(nbrFaceCells, facei)
238  {
239  valueDiff -= p[nbrFaceCells[facei]];
240  }
241 
242  Info<< "Pressure difference = " << valueDiff << endl;
243  }
244 
245  if ((mag(valueDiff) > mag(minThresholdValue_)) || baffleActivated_)
246  {
247  openFraction_ =
248  max(
249  min(
250  openFraction_
251  + min
252  (
253  this->db().time().deltaT().value()/openingTime_,
254  maxOpenFractionDelta_
255  )*(orientation_),
256  1 - 1e-6
257  ),
258  1e-6
259  );
260 
261  baffleActivated_ = true;
262  }
263  else
264  {
265  openFraction_ = max(min(1 - 1e-6, openFraction_), 1e-6);
266  }
267 
268  Info<< "Open fraction = " << openFraction_ << endl;
269 
270  vectorField::subField Sfw = patch().patch().faceAreas();
271  vectorField newSfw((1 - openFraction_)*initWallSf_);
272  forAll(Sfw, facei)
273  {
274  Sfw[facei] = newSfw[facei];
275  }
276  const_cast<scalarField&>(patch().magSf()) = mag(patch().Sf());
277 
278  // Update owner side of cyclic
279  const_cast<vectorField&>(cyclicPatch.Sf()) =
280  openFraction_*initCyclicSf_;
281 
282  const_cast<scalarField&>(cyclicPatch.magSf()) =
283  mag(cyclicPatch.Sf());
284 
285  // Update neighbour side of cyclic
286  const_cast<vectorField&>(nbrPatch.Sf()) =
287  openFraction_*nbrCyclicSf_;
288 
289  const_cast<scalarField&>(nbrPatch.magSf()) =
290  mag(nbrPatch.Sf());
291 
292  curTimeIndex_ = this->db().time().timeIndex();
293  }
294 
295  fixedValueFvPatchVectorField::updateCoeffs();
296 }
297 
298 
300 write(Ostream& os) const
301 {
303  writeEntryIfDifferent<word>(os, "p", "p", pName_);
304  writeEntry(os, "cyclicPatch", cyclicPatchName_);
305  writeEntry(os, "orientation", orientation_);
306  writeEntry(os, "openingTime", openingTime_);
307  writeEntry(os, "maxOpenFractionDelta", maxOpenFractionDelta_);
308  writeEntry(os, "openFraction", openFraction_);
309  writeEntry(os, "minThresholdValue", minThresholdValue_);
310  writeEntry(os, "forceBased", fBased_);
311  writeEntry(os, "value", *this);
312 }
313 
314 
315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316 
317 namespace Foam
318 {
320  (
323  );
324 }
325 
326 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Pre-declare related SubField type.
Definition: SubField.H:63
This boundary condition is applied to the flow velocity, to simulate the opening or closure of a baff...
activePressureForceBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchVectorField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:120
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:142
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:162
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:136
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:313
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const unitConversion unitFraction
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict
volScalarField & p
Foam::surfaceFields.