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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "cyclicFvPatch.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  fixedValueFvPatchVectorField(p, iF),
44  pName_("p"),
45  cyclicPatchName_(),
46  cyclicPatchLabel_(-1),
47  orientation_(1),
48  initWallSf_(0),
49  initCyclicSf_(0),
50  nbrCyclicSf_(0),
51  openFraction_(0),
52  openingTime_(0),
53  maxOpenFractionDelta_(0),
54  curTimeIndex_(-1)
55 {}
56 
57 
60 (
62  const fvPatch& p,
64  const fvPatchFieldMapper& mapper
65 )
66 :
67  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
68  pName_(ptf.pName_),
69  cyclicPatchName_(ptf.cyclicPatchName_),
70  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
71  orientation_(ptf.orientation_),
72  initWallSf_(ptf.initWallSf_),
73  initCyclicSf_(ptf.initCyclicSf_),
74  nbrCyclicSf_(ptf.nbrCyclicSf_),
75  openFraction_(ptf.openFraction_),
76  openingTime_(ptf.openingTime_),
77  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
78  curTimeIndex_(-1)
79 {}
80 
81 
84 (
85  const fvPatch& p,
87  const dictionary& dict
88 )
89 :
90  fixedValueFvPatchVectorField(p, iF, dict, false),
91  pName_(dict.lookupOrDefault<word>("p", "p")),
92  cyclicPatchName_(dict.lookup("cyclicPatch")),
93  cyclicPatchLabel_(p.patch().boundaryMesh().findPatchID(cyclicPatchName_)),
94  orientation_(dict.get<label>("orientation")),
95  initWallSf_(p.Sf()),
96  initCyclicSf_(p.boundaryMesh()[cyclicPatchLabel_].Sf()),
97  nbrCyclicSf_
98  (
99  refCast<const cyclicFvPatch>
100  (
101  p.boundaryMesh()[cyclicPatchLabel_],
102  dict
103  ).neighbFvPatch().Sf()
104  ),
105  openFraction_(dict.get<scalar>("openFraction")),
106  openingTime_(dict.get<scalar>("openingTime")),
107  maxOpenFractionDelta_(dict.get<scalar>("maxOpenFractionDelta")),
108  curTimeIndex_(-1)
109 {
110  fvPatchVectorField::operator=(Zero);
111 }
112 
113 
116 (
118 )
119 :
120  fixedValueFvPatchVectorField(ptf),
121  pName_(ptf.pName_),
122  cyclicPatchName_(ptf.cyclicPatchName_),
123  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
124  orientation_(ptf.orientation_),
125  initWallSf_(ptf.initWallSf_),
126  initCyclicSf_(ptf.initCyclicSf_),
127  nbrCyclicSf_(ptf.nbrCyclicSf_),
128  openFraction_(ptf.openFraction_),
129  openingTime_(ptf.openingTime_),
130  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
131  curTimeIndex_(-1)
132 {}
133 
134 
137 (
140 )
141 :
142  fixedValueFvPatchVectorField(ptf, iF),
143  pName_(ptf.pName_),
144  cyclicPatchName_(ptf.cyclicPatchName_),
145  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
146  orientation_(ptf.orientation_),
147  initWallSf_(ptf.initWallSf_),
148  initCyclicSf_(ptf.initCyclicSf_),
149  nbrCyclicSf_(ptf.nbrCyclicSf_),
150  openFraction_(ptf.openFraction_),
151  openingTime_(ptf.openingTime_),
152  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
153  curTimeIndex_(-1)
154 {}
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
160 (
161  const fvPatchFieldMapper& m
162 )
163 {
164  fixedValueFvPatchVectorField::autoMap(m);
165 
166  //- Note: cannot map field from cyclic patch anyway so just recalculate
167  // Areas should be consistent when doing autoMap except in case of
168  // topo changes.
169  //- Note: we don't want to use Sf here since triggers rebuilding of
170  // fvMesh::S() which will give problems when mapped (since already
171  // on new mesh)
172  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
173  initWallSf_ = patch().patchSlice(areas);
174  initCyclicSf_ = patch().boundaryMesh()
175  [
176  cyclicPatchLabel_
177  ].patchSlice(areas);
178  nbrCyclicSf_ = refCast<const cyclicFvPatch>
179  (
180  patch().boundaryMesh()
181  [
182  cyclicPatchLabel_
183  ]
184  ).neighbFvPatch().patch().patchSlice(areas);
185 }
186 
187 
189 (
190  const fvPatchVectorField& ptf,
191  const labelList& addr
192 )
193 {
194  fixedValueFvPatchVectorField::rmap(ptf, addr);
195 
196  // See autoMap.
197  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
198  initWallSf_ = patch().patchSlice(areas);
199  initCyclicSf_ = patch().boundaryMesh()
200  [
201  cyclicPatchLabel_
202  ].patchSlice(areas);
203  nbrCyclicSf_ = refCast<const cyclicFvPatch>
204  (
205  patch().boundaryMesh()
206  [
207  cyclicPatchLabel_
208  ]
209  ).neighbFvPatch().patch().patchSlice(areas);
210 }
211 
212 
214 {
215  if (updated())
216  {
217  return;
218  }
219 
220  // Execute the change to the openFraction only once per time-step
221  if (curTimeIndex_ != this->db().time().timeIndex())
222  {
223  const volScalarField& p = db().lookupObject<volScalarField>
224  (
225  pName_
226  );
227 
228  const fvPatch& cyclicPatch = patch().boundaryMesh()[cyclicPatchLabel_];
229  const labelList& cyclicFaceCells = cyclicPatch.patch().faceCells();
230  const fvPatch& nbrPatch = refCast<const cyclicFvPatch>
231  (
232  cyclicPatch
233  ).neighbFvPatch();
234  const labelList& nbrFaceCells = nbrPatch.patch().faceCells();
235 
236  scalar forceDiff = 0;
237 
238  // Add this side
239  forAll(cyclicFaceCells, facei)
240  {
241  forceDiff += p[cyclicFaceCells[facei]]*mag(initCyclicSf_[facei]);
242  }
243 
244  // Remove other side
245  forAll(nbrFaceCells, facei)
246  {
247  forceDiff -= p[nbrFaceCells[facei]]*mag(nbrCyclicSf_[facei]);
248  }
249 
250  openFraction_ =
251  max
252  (
253  min
254  (
255  openFraction_
256  + min
257  (
258  this->db().time().deltaTValue()/openingTime_,
259  maxOpenFractionDelta_
260  )
261  *(orientation_*sign(forceDiff)),
262  1 - 1e-6
263  ),
264  1e-6
265  );
266 
267  Info<< "openFraction = " << openFraction_ << endl;
268 
269  vectorField::subField Sfw = this->patch().patch().faceAreas();
270  const vectorField newSfw((1 - openFraction_)*initWallSf_);
271  forAll(Sfw, facei)
272  {
273  Sfw[facei] = newSfw[facei];
274  }
275  const_cast<scalarField&>(patch().magSf()) = mag(patch().Sf());
276 
277  // Update owner side of cyclic
278  const_cast<vectorField&>(cyclicPatch.Sf()) =
279  openFraction_*initCyclicSf_;
280  const_cast<scalarField&>(cyclicPatch.magSf()) =
281  mag(cyclicPatch.Sf());
282  // Update neighbour side of cyclic
283  const_cast<vectorField&>(nbrPatch.Sf()) =
284  openFraction_*nbrCyclicSf_;
285  const_cast<scalarField&>(nbrPatch.magSf()) =
286  mag(nbrPatch.Sf());
287 
288  curTimeIndex_ = this->db().time().timeIndex();
289  }
290 
291  fixedValueFvPatchVectorField::updateCoeffs();
292 }
293 
294 
296 {
298  os.writeEntryIfDifferent<word>("p", "p", pName_);
299  os.writeEntry("cyclicPatch", cyclicPatchName_);
300  os.writeEntry("orientation", orientation_);
301  os.writeEntry("openingTime", openingTime_);
302  os.writeEntry("maxOpenFractionDelta", maxOpenFractionDelta_);
303  os.writeEntry("openFraction", openFraction_);
304  writeEntry("value", os);
305 }
306 
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 namespace Foam
311 {
313  (
316  );
317 }
318 
319 
320 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::fvPatch::Sf
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:126
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:231
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:62
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::activeBaffleVelocityFvPatchVectorField::rmap
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: activeBaffleVelocityFvPatchVectorField.C:189
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
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:295
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::SubField
SubField is a Field obtained as a section of another Field.
Definition: Field.H:64
Foam::label
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:62
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::activeBaffleVelocityFvPatchVectorField::activeBaffleVelocityFvPatchVectorField
activeBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: activeBaffleVelocityFvPatchVectorField.C:38
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::activeBaffleVelocityFvPatchVectorField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: activeBaffleVelocityFvPatchVectorField.C:160
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
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
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:121
Foam::activeBaffleVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: activeBaffleVelocityFvPatchVectorField.C:213
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:361
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:219
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:145
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:61
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
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:132
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