activeBaffleVelocityFvPatchVectorField.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-2021 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 {}
54 
55 
58 (
59  const fvPatch& p,
61  const dictionary& dict
62 )
63 :
64  fixedValueFvPatchVectorField(p, iF, dict, false),
65  pName_(dict.lookupOrDefault<word>("p", "p")),
66  cyclicPatchName_(dict.lookup("cyclicPatch")),
67  cyclicPatchLabel_(p.patch().boundaryMesh().findPatchID(cyclicPatchName_)),
68  orientation_(dict.lookup<label>("orientation")),
69  initWallSf_(p.Sf()),
70  initCyclicSf_(p.boundaryMesh()[cyclicPatchLabel_].Sf()),
71  nbrCyclicSf_
72  (
73  refCast<const cyclicFvPatch>
74  (
75  p.boundaryMesh()[cyclicPatchLabel_]
76  ).neighbFvPatch().Sf()
77  ),
78  openFraction_(dict.lookup<scalar>("openFraction")),
79  openingTime_(dict.lookup<scalar>("openingTime")),
80  maxOpenFractionDelta_(dict.lookup<scalar>("maxOpenFractionDelta")),
81  curTimeIndex_(-1)
82 {
83  fvPatchVectorField::operator=(Zero);
84 }
85 
86 
89 (
91  const fvPatch& p,
93  const fvPatchFieldMapper& mapper
94 )
95 :
96  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
97  pName_(ptf.pName_),
98  cyclicPatchName_(ptf.cyclicPatchName_),
99  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
100  orientation_(ptf.orientation_),
101  initWallSf_(ptf.initWallSf_),
102  initCyclicSf_(ptf.initCyclicSf_),
103  nbrCyclicSf_(ptf.nbrCyclicSf_),
104  openFraction_(ptf.openFraction_),
105  openingTime_(ptf.openingTime_),
106  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
107  curTimeIndex_(-1)
108 {}
109 
110 
113 (
116 )
117 :
118  fixedValueFvPatchVectorField(ptf, iF),
119  pName_(ptf.pName_),
120  cyclicPatchName_(ptf.cyclicPatchName_),
121  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
122  orientation_(ptf.orientation_),
123  initWallSf_(ptf.initWallSf_),
124  initCyclicSf_(ptf.initCyclicSf_),
125  nbrCyclicSf_(ptf.nbrCyclicSf_),
126  openFraction_(ptf.openFraction_),
127  openingTime_(ptf.openingTime_),
128  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
129  curTimeIndex_(-1)
130 {}
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
136 (
137  const fvPatchFieldMapper& m
138 )
139 {
140  fixedValueFvPatchVectorField::autoMap(m);
141 
142  //- Note: cannot map field from cyclic patch anyway so just recalculate
143  // Areas should be consistent when doing autoMap except in case of
144  // topo 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  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
149  initWallSf_ = patch().patchSlice(areas);
150  initCyclicSf_ = patch().boundaryMesh()
151  [
152  cyclicPatchLabel_
153  ].patchSlice(areas);
154  nbrCyclicSf_ = refCast<const cyclicFvPatch>
155  (
156  patch().boundaryMesh()
157  [
158  cyclicPatchLabel_
159  ]
160  ).neighbFvPatch().patch().patchSlice(areas);
161 }
162 
163 
165 (
166  const fvPatchVectorField& ptf,
167  const labelList& addr
168 )
169 {
170  fixedValueFvPatchVectorField::rmap(ptf, addr);
171 
172  // See autoMap.
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 
196  // Execute the change to the openFraction only once per time-step
197  if (curTimeIndex_ != this->db().time().timeIndex())
198  {
199  const volScalarField& p = db().lookupObject<volScalarField>
200  (
201  pName_
202  );
203 
204  const fvPatch& cyclicPatch = patch().boundaryMesh()[cyclicPatchLabel_];
205  const labelList& cyclicFaceCells = cyclicPatch.patch().faceCells();
206  const fvPatch& nbrPatch = refCast<const cyclicFvPatch>
207  (
208  cyclicPatch
209  ).neighbFvPatch();
210  const labelList& nbrFaceCells = nbrPatch.patch().faceCells();
211 
212  scalar forceDiff = 0;
213 
214  // Add this side
215  forAll(cyclicFaceCells, facei)
216  {
217  forceDiff += p[cyclicFaceCells[facei]]*mag(initCyclicSf_[facei]);
218  }
219 
220  // Remove other side
221  forAll(nbrFaceCells, facei)
222  {
223  forceDiff -= p[nbrFaceCells[facei]]*mag(nbrCyclicSf_[facei]);
224  }
225 
226  openFraction_ =
227  max
228  (
229  min
230  (
231  openFraction_
232  + min
233  (
234  this->db().time().deltaTValue()/openingTime_,
235  maxOpenFractionDelta_
236  )
237  *(orientation_*sign(forceDiff)),
238  1 - 1e-6
239  ),
240  1e-6
241  );
242 
243  Info<< "openFraction = " << openFraction_ << endl;
244 
245  vectorField::subField Sfw = this->patch().patch().faceAreas();
246  const vectorField newSfw((1 - openFraction_)*initWallSf_);
247  forAll(Sfw, facei)
248  {
249  Sfw[facei] = newSfw[facei];
250  }
251  const_cast<scalarField&>(patch().magSf()) = mag(patch().Sf());
252 
253  // Update owner side of cyclic
254  const_cast<vectorField&>(cyclicPatch.Sf()) =
255  openFraction_*initCyclicSf_;
256  const_cast<scalarField&>(cyclicPatch.magSf()) =
257  mag(cyclicPatch.Sf());
258  // Update neighbour side of cyclic
259  const_cast<vectorField&>(nbrPatch.Sf()) =
260  openFraction_*nbrCyclicSf_;
261  const_cast<scalarField&>(nbrPatch.magSf()) =
262  mag(nbrPatch.Sf());
263 
264  curTimeIndex_ = this->db().time().timeIndex();
265  }
266 
267  fixedValueFvPatchVectorField::updateCoeffs();
268 }
269 
270 
272 {
274  writeEntryIfDifferent<word>(os, "p", "p", pName_);
275  writeEntry(os, "cyclicPatch", cyclicPatchName_);
276  writeEntry(os, "orientation", orientation_);
277  writeEntry(os, "openingTime", openingTime_);
278  writeEntry(os, "maxOpenFractionDelta", maxOpenFractionDelta_);
279  writeEntry(os, "openFraction", openFraction_);
280  writeEntryIfDifferent<word>(os, "p", "p", pName_);
281  writeEntry(os, "value", *this);
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 namespace Foam
288 {
290  (
293  );
294 }
295 
296 
297 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:180
dimensionedScalar sign(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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 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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:291
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:62
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.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:260
Pre-declare related SubField type.
Definition: Field.H:60
Macros for easy insertion into run-time selection tables.
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:334
A class for handling words, derived from string.
Definition: word.H:59
activeBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:138
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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
This velocity boundary condition simulates the opening of a baffle due to local flow conditions...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844