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-2023 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().findPatchID(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")),
57  openingTime_(dict.lookup<scalar>("openingTime")),
58  maxOpenFractionDelta_(dict.lookup<scalar>("maxOpenFractionDelta")),
59  curTimeIndex_(-1)
60 {
62 }
63 
64 
67 (
69  const fvPatch& p,
71  const fvPatchFieldMapper& mapper
72 )
73 :
74  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
75  pName_(ptf.pName_),
76  cyclicPatchName_(ptf.cyclicPatchName_),
77  cyclicPatchLabel_(ptf.cyclicPatchLabel_),
78  orientation_(ptf.orientation_),
79  initWallSf_(ptf.initWallSf_),
80  initCyclicSf_(ptf.initCyclicSf_),
81  nbrCyclicSf_(ptf.nbrCyclicSf_),
82  openFraction_(ptf.openFraction_),
83  openingTime_(ptf.openingTime_),
84  maxOpenFractionDelta_(ptf.maxOpenFractionDelta_),
85  curTimeIndex_(-1)
86 {}
87 
88 
91 (
94 )
95 :
96  fixedValueFvPatchVectorField(ptf, iF),
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 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 (
115  const fvPatchVectorField& ptf,
116  const fvPatchFieldMapper& mapper
117 )
118 {
119  fixedValueFvPatchVectorField::map(ptf, mapper);
120 
121  //- Note: cannot map field from cyclic patch anyway so just recalculate
122  // Areas should be consistent when doing map except in case of topo
123  // changes.
124  //- Note: we don't want to use Sf here since triggers rebuilding of
125  // fvMesh::S() which will give problems when mapped (since already
126  // on new mesh)
127 
128  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
129  initWallSf_ = patch().patchSlice(areas);
130  initCyclicSf_ = patch().boundaryMesh()
131  [
132  cyclicPatchLabel_
133  ].patchSlice(areas);
134  nbrCyclicSf_ = refCast<const cyclicFvPatch>
135  (
136  patch().boundaryMesh()
137  [
138  cyclicPatchLabel_
139  ]
140  ).neighbFvPatch().patch().patchSlice(areas);
141 }
142 
143 
145 (
146  const fvPatchVectorField& ptf
147 )
148 {
149  fixedValueFvPatchVectorField::reset(ptf);
150 
151  // See rmap.
152  const vectorField& areas = patch().boundaryMesh().mesh().faceAreas();
153  initWallSf_ = patch().patchSlice(areas);
154  initCyclicSf_ = patch().boundaryMesh()
155  [
156  cyclicPatchLabel_
157  ].patchSlice(areas);
158  nbrCyclicSf_ = refCast<const cyclicFvPatch>
159  (
160  patch().boundaryMesh()
161  [
162  cyclicPatchLabel_
163  ]
164  ).neighbFvPatch().patch().patchSlice(areas);
165 }
166 
167 
169 {
170  if (updated())
171  {
172  return;
173  }
174 
175  // Execute the change to the openFraction only once per time-step
176  if (curTimeIndex_ != this->db().time().timeIndex())
177  {
178  const volScalarField& p = db().lookupObject<volScalarField>
179  (
180  pName_
181  );
182 
183  const fvPatch& cyclicPatch = patch().boundaryMesh()[cyclicPatchLabel_];
184  const labelList& cyclicFaceCells = cyclicPatch.patch().faceCells();
185  const fvPatch& nbrPatch = refCast<const cyclicFvPatch>
186  (
187  cyclicPatch
188  ).neighbFvPatch();
189  const labelList& nbrFaceCells = nbrPatch.patch().faceCells();
190 
191  scalar forceDiff = 0;
192 
193  // Add this side
194  forAll(cyclicFaceCells, facei)
195  {
196  forceDiff += p[cyclicFaceCells[facei]]*mag(initCyclicSf_[facei]);
197  }
198 
199  // Remove other side
200  forAll(nbrFaceCells, facei)
201  {
202  forceDiff -= p[nbrFaceCells[facei]]*mag(nbrCyclicSf_[facei]);
203  }
204 
205  openFraction_ =
206  max
207  (
208  min
209  (
210  openFraction_
211  + min
212  (
213  this->db().time().deltaTValue()/openingTime_,
214  maxOpenFractionDelta_
215  )
216  *(orientation_*sign(forceDiff)),
217  1 - 1e-6
218  ),
219  1e-6
220  );
221 
222  Info<< "openFraction = " << openFraction_ << endl;
223 
224  vectorField::subField Sfw = this->patch().patch().faceAreas();
225  const vectorField newSfw((1 - openFraction_)*initWallSf_);
226  forAll(Sfw, facei)
227  {
228  Sfw[facei] = newSfw[facei];
229  }
230  const_cast<scalarField&>(patch().magSf()) = mag(patch().Sf());
231 
232  // Update owner side of cyclic
233  const_cast<vectorField&>(cyclicPatch.Sf()) =
234  openFraction_*initCyclicSf_;
235  const_cast<scalarField&>(cyclicPatch.magSf()) =
236  mag(cyclicPatch.Sf());
237  // Update neighbour side of cyclic
238  const_cast<vectorField&>(nbrPatch.Sf()) =
239  openFraction_*nbrCyclicSf_;
240  const_cast<scalarField&>(nbrPatch.magSf()) =
241  mag(nbrPatch.Sf());
242 
243  curTimeIndex_ = this->db().time().timeIndex();
244  }
245 
246  fixedValueFvPatchVectorField::updateCoeffs();
247 }
248 
249 
251 {
253  writeEntryIfDifferent<word>(os, "p", "p", pName_);
254  writeEntry(os, "cyclicPatch", cyclicPatchName_);
255  writeEntry(os, "orientation", orientation_);
256  writeEntry(os, "openingTime", openingTime_);
257  writeEntry(os, "maxOpenFractionDelta", maxOpenFractionDelta_);
258  writeEntry(os, "openFraction", openFraction_);
259  writeEntryIfDifferent<word>(os, "p", "p", pName_);
260  writeEntry(os, "value", *this);
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
266 namespace Foam
267 {
269  (
272  );
273 }
274 
275 
276 // ************************************************************************* //
#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 map(const fvPatchVectorField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given 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:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:139
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:142
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:181
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:105
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:111
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict
volScalarField & p
Foam::surfaceFields.