activeBaffleVelocityFvPatchVectorField.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 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  orientation_(1),
49  initWallSf_(0),
50  initCyclicSf_(0),
51  nbrCyclicSf_(0),
52  openFraction_(0),
53  openingTime_(0),
54  maxOpenFractionDelta_(0),
55  curTimeIndex_(-1)
56 {}
57 
58 
61 (
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
68  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
69  pName_(ptf.pName_),
70  cyclicPatchName_(ptf.cyclicPatchName_),
71  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
72  orientation_(ptf.orientation_),
73  initWallSf_(ptf.initWallSf_),
74  initCyclicSf_(ptf.initCyclicSf_),
75  nbrCyclicSf_(ptf.nbrCyclicSf_),
76  openFraction_(ptf.openFraction_),
77  openingTime_(ptf.openingTime_),
78  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
79  curTimeIndex_(-1)
80 {}
81 
82 
85 (
86  const fvPatch& p,
88  const dictionary& dict
89 )
90 :
91  fixedValueFvPatchVectorField(p, iF, dict, false),
92  pName_(dict.getOrDefault<word>("p", "p")),
93  cyclicPatchName_(dict.lookup("cyclicPatch")),
94  cyclicPatchLabel_(p.patch().boundaryMesh().findPatchID(cyclicPatchName_)),
95  orientation_(dict.get<label>("orientation")),
96  initWallSf_(p.Sf()),
97  initCyclicSf_(p.boundaryMesh()[cyclicPatchLabel_].Sf()),
98  nbrCyclicSf_
99  (
100  refCast<const cyclicFvPatch>
101  (
102  p.boundaryMesh()[cyclicPatchLabel_],
103  dict
104  ).neighbFvPatch().Sf()
105  ),
106  openFraction_(dict.get<scalar>("openFraction")),
107  openingTime_(dict.get<scalar>("openingTime")),
108  maxOpenFractionDelta_(dict.get<scalar>("maxOpenFractionDelta")),
109  curTimeIndex_(-1)
110 {
111  fvPatchVectorField::operator=(Zero);
112 }
113 
114 
117 (
119 )
120 :
121  fixedValueFvPatchVectorField(ptf),
122  pName_(ptf.pName_),
123  cyclicPatchName_(ptf.cyclicPatchName_),
124  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
125  orientation_(ptf.orientation_),
126  initWallSf_(ptf.initWallSf_),
127  initCyclicSf_(ptf.initCyclicSf_),
128  nbrCyclicSf_(ptf.nbrCyclicSf_),
129  openFraction_(ptf.openFraction_),
130  openingTime_(ptf.openingTime_),
131  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
132  curTimeIndex_(-1)
133 {}
134 
135 
138 (
141 )
142 :
143  fixedValueFvPatchVectorField(ptf, iF),
144  pName_(ptf.pName_),
145  cyclicPatchName_(ptf.cyclicPatchName_),
146  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
147  orientation_(ptf.orientation_),
148  initWallSf_(ptf.initWallSf_),
149  initCyclicSf_(ptf.initCyclicSf_),
150  nbrCyclicSf_(ptf.nbrCyclicSf_),
151  openFraction_(ptf.openFraction_),
152  openingTime_(ptf.openingTime_),
153  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
154  curTimeIndex_(-1)
155 {}
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
161 (
162  const fvPatchFieldMapper& m
163 )
164 {
165  fixedValueFvPatchVectorField::autoMap(m);
166 
167  //- Note: cannot map field from cyclic patch anyway so just recalculate
168  // Areas should be consistent when doing autoMap except in case of
169  // topo changes.
170  //- Note: we don't want to use Sf here since triggers rebuilding of
171  // fvMesh::S() which will give problems when mapped (since already
172  // on new mesh)
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  const fvPatchVectorField& ptf,
192  const labelList& addr
193 )
194 {
195  fixedValueFvPatchVectorField::rmap(ptf, addr);
196 
197  // See autoMap.
198  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
199  initWallSf_ = patch().patchSlice(areas);
200  initCyclicSf_ = patch().boundaryMesh()
201  [
202  cyclicPatchLabel_
203  ].patchSlice(areas);
204  nbrCyclicSf_ = refCast<const cyclicFvPatch>
205  (
206  patch().boundaryMesh()
207  [
208  cyclicPatchLabel_
209  ]
210  ).neighbFvPatch().patch().patchSlice(areas);
211 }
212 
213 
215 {
216  if (updated())
217  {
218  return;
219  }
220 
221  // Execute the change to the openFraction only once per time-step
222  if (curTimeIndex_ != this->db().time().timeIndex())
223  {
224  const volScalarField& p = db().lookupObject<volScalarField>
225  (
226  pName_
227  );
228 
229  const fvPatch& cyclicPatch = patch().boundaryMesh()[cyclicPatchLabel_];
230  const labelList& cyclicFaceCells = cyclicPatch.patch().faceCells();
231  const fvPatch& nbrPatch = refCast<const cyclicFvPatch>
232  (
233  cyclicPatch
234  ).neighbFvPatch();
235  const labelList& nbrFaceCells = nbrPatch.patch().faceCells();
236 
237  scalar forceDiff = 0;
238 
239  // Add this side
240  forAll(cyclicFaceCells, facei)
241  {
242  forceDiff += p[cyclicFaceCells[facei]]*mag(initCyclicSf_[facei]);
243  }
244 
245  // Remove other side
246  forAll(nbrFaceCells, facei)
247  {
248  forceDiff -= p[nbrFaceCells[facei]]*mag(nbrCyclicSf_[facei]);
249  }
250 
251  openFraction_ =
252  max
253  (
254  min
255  (
256  openFraction_
257  + min
258  (
259  this->db().time().deltaTValue()/openingTime_,
260  maxOpenFractionDelta_
261  )
262  *(orientation_*sign(forceDiff)),
263  1 - 1e-6
264  ),
265  1e-6
266  );
267 
268  Info<< "openFraction = " << openFraction_ << endl;
269 
270  vectorField::subField Sfw = this->patch().patch().faceAreas();
271  const 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  const_cast<scalarField&>(cyclicPatch.magSf()) =
282  mag(cyclicPatch.Sf());
283  // Update neighbour side of cyclic
284  const_cast<vectorField&>(nbrPatch.Sf()) =
285  openFraction_*nbrCyclicSf_;
286  const_cast<scalarField&>(nbrPatch.magSf()) =
287  mag(nbrPatch.Sf());
288 
289  curTimeIndex_ = this->db().time().timeIndex();
290  }
291 
292  fixedValueFvPatchVectorField::updateCoeffs();
293 }
294 
295 
297 {
299  os.writeEntryIfDifferent<word>("p", "p", pName_);
300  os.writeEntry("cyclicPatch", cyclicPatchName_);
301  os.writeEntry("orientation", orientation_);
302  os.writeEntry("openingTime", openingTime_);
303  os.writeEntry("maxOpenFractionDelta", maxOpenFractionDelta_);
304  os.writeEntry("openFraction", openFraction_);
305  writeEntry("value", os);
306 }
307 
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 namespace Foam
312 {
314  (
317  );
318 }
319 
320 
321 // ************************************************************************* //
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
activeBaffleVelocityFvPatchVectorField.H
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::activeBaffleVelocityFvPatchVectorField::rmap
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: activeBaffleVelocityFvPatchVectorField.C:190
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
surfaceFields.H
Foam::surfaceFields.
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
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::activeBaffleVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: activeBaffleVelocityFvPatchVectorField.C:296
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::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::activeBaffleVelocityFvPatchVectorField::activeBaffleVelocityFvPatchVectorField
activeBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: activeBaffleVelocityFvPatchVectorField.C:39
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::activeBaffleVelocityFvPatchVectorField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: activeBaffleVelocityFvPatchVectorField.C:161
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::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::activeBaffleVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: activeBaffleVelocityFvPatchVectorField.C:214
os
OBJstream os(runTime.globalPath()/outputName)
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::activeBaffleVelocityFvPatchVectorField
This velocity boundary condition simulates the opening of a baffle due to local flow conditions,...
Definition: activeBaffleVelocityFvPatchVectorField.H:162
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::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
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::GeometricField< scalar, fvPatchField, volMesh >
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54