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-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_(p.Sf()),
48  initCyclicSf_(p.boundaryMesh()[cyclicPatchLabel_].Sf()),
49  nbrCyclicSf_
50  (
51  refCast<const cyclicFvPatch>
52  (
53  p.boundaryMesh()[cyclicPatchLabel_]
54  ).neighbFvPatch().Sf()
55  ),
56  openFraction_(dict.lookup<scalar>("openFraction", unitFraction)),
57  openingTime_(dict.lookup<scalar>("openingTime", dimTime)),
58  maxOpenFractionDelta_
59  (
60  dict.lookup<scalar>("maxOpenFractionDelta", unitFraction)
61  ),
62  curTimeIndex_(-1)
63 {
65 }
66 
67 
70 (
72  const fvPatch& p,
74  const fieldMapper& mapper
75 )
76 :
77  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
78  pName_(ptf.pName_),
79  cyclicPatchName_(ptf.cyclicPatchName_),
80  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
81  orientation_(ptf.orientation_),
82  initWallSf_(ptf.initWallSf_),
83  initCyclicSf_(ptf.initCyclicSf_),
84  nbrCyclicSf_(ptf.nbrCyclicSf_),
85  openFraction_(ptf.openFraction_),
86  openingTime_(ptf.openingTime_),
87  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
88  curTimeIndex_(-1)
89 {}
90 
91 
94 (
97 )
98 :
99  fixedValueFvPatchVectorField(ptf, iF),
100  pName_(ptf.pName_),
101  cyclicPatchName_(ptf.cyclicPatchName_),
102  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
103  orientation_(ptf.orientation_),
104  initWallSf_(ptf.initWallSf_),
105  initCyclicSf_(ptf.initCyclicSf_),
106  nbrCyclicSf_(ptf.nbrCyclicSf_),
107  openFraction_(ptf.openFraction_),
108  openingTime_(ptf.openingTime_),
109  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
110  curTimeIndex_(-1)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 (
118  const fvPatchVectorField& ptf,
119  const fieldMapper& mapper
120 )
121 {
122  fixedValueFvPatchVectorField::map(ptf, mapper);
123 
124  //- Note: cannot map field from cyclic patch anyway so just recalculate
125  // Areas should be consistent when doing map except in case of topo
126  // changes.
127  //- Note: we don't want to use Sf here since triggers rebuilding of
128  // fvMesh::S() which will give problems when mapped (since already
129  // on new mesh)
130 
131  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
132  initWallSf_ = patch().patchSlice(areas);
133  initCyclicSf_ = patch().boundaryMesh()
134  [
135  cyclicPatchLabel_
136  ].patchSlice(areas);
137  nbrCyclicSf_ = refCast<const cyclicFvPatch>
138  (
139  patch().boundaryMesh()
140  [
141  cyclicPatchLabel_
142  ]
143  ).neighbFvPatch().patch().patchSlice(areas);
144 }
145 
146 
148 (
149  const fvPatchVectorField& ptf
150 )
151 {
152  fixedValueFvPatchVectorField::reset(ptf);
153 
154  // See rmap.
155  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
156  initWallSf_ = patch().patchSlice(areas);
157  initCyclicSf_ = patch().boundaryMesh()
158  [
159  cyclicPatchLabel_
160  ].patchSlice(areas);
161  nbrCyclicSf_ = refCast<const cyclicFvPatch>
162  (
163  patch().boundaryMesh()
164  [
165  cyclicPatchLabel_
166  ]
167  ).neighbFvPatch().patch().patchSlice(areas);
168 }
169 
170 
172 {
173  if (updated())
174  {
175  return;
176  }
177 
178  // Execute the change to the openFraction only once per time-step
179  if (curTimeIndex_ != this->db().time().timeIndex())
180  {
181  const volScalarField& p = db().lookupObject<volScalarField>
182  (
183  pName_
184  );
185 
186  const fvPatch& cyclicPatch = patch().boundaryMesh()[cyclicPatchLabel_];
187  const labelList& cyclicFaceCells = cyclicPatch.patch().faceCells();
188  const fvPatch& nbrPatch = refCast<const cyclicFvPatch>
189  (
190  cyclicPatch
191  ).neighbFvPatch();
192  const labelList& nbrFaceCells = nbrPatch.patch().faceCells();
193 
194  scalar forceDiff = 0;
195 
196  // Add this side
197  forAll(cyclicFaceCells, facei)
198  {
199  forceDiff += p[cyclicFaceCells[facei]]*mag(initCyclicSf_[facei]);
200  }
201 
202  // Remove other side
203  forAll(nbrFaceCells, facei)
204  {
205  forceDiff -= p[nbrFaceCells[facei]]*mag(nbrCyclicSf_[facei]);
206  }
207 
208  openFraction_ =
209  max
210  (
211  min
212  (
213  openFraction_
214  + min
215  (
216  this->db().time().deltaTValue()/openingTime_,
217  maxOpenFractionDelta_
218  )
219  *(orientation_*sign(forceDiff)),
220  1 - 1e-6
221  ),
222  1e-6
223  );
224 
225  Info<< "openFraction = " << openFraction_ << endl;
226 
227  vectorField::subField Sfw = this->patch().patch().faceAreas();
228  const vectorField newSfw((1 - openFraction_)*initWallSf_);
229  forAll(Sfw, facei)
230  {
231  Sfw[facei] = newSfw[facei];
232  }
233  const_cast<scalarField&>(patch().magSf()) = mag(patch().Sf());
234 
235  // Update owner side of cyclic
236  const_cast<vectorField&>(cyclicPatch.Sf()) =
237  openFraction_*initCyclicSf_;
238  const_cast<scalarField&>(cyclicPatch.magSf()) =
239  mag(cyclicPatch.Sf());
240  // Update neighbour side of cyclic
241  const_cast<vectorField&>(nbrPatch.Sf()) =
242  openFraction_*nbrCyclicSf_;
243  const_cast<scalarField&>(nbrPatch.magSf()) =
244  mag(nbrPatch.Sf());
245 
246  curTimeIndex_ = this->db().time().timeIndex();
247  }
248 
249  fixedValueFvPatchVectorField::updateCoeffs();
250 }
251 
252 
254 {
256  writeEntryIfDifferent<word>(os, "p", "p", pName_);
257  writeEntry(os, "cyclicPatch", cyclicPatchName_);
258  writeEntry(os, "orientation", orientation_);
259  writeEntry(os, "openingTime", openingTime_);
260  writeEntry(os, "maxOpenFractionDelta", maxOpenFractionDelta_);
261  writeEntry(os, "openFraction", openFraction_);
262  writeEntryIfDifferent<word>(os, "p", "p", pName_);
263  writeEntry(os, "value", *this);
264 }
265 
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 namespace Foam
270 {
272  (
275  );
276 }
277 
278 
279 // ************************************************************************* //
#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 velocity boundary condition simulates the opening of a baffle due to local flow conditions,...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
activeBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
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.
Cyclic-plane patch.
Definition: cyclicFvPatch.H:55
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
dimensionedScalar sign(const dimensionedScalar &ds)
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
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:129
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.