wallCellWallFunctionFvPatchScalarField.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 
27 #include "fvMatrix.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 Foam::scalar Foam::wallCellWallFunctionFvPatchScalarField::tol_ = 1e-1;
33 
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
38 {
39  volScalarField& vf =
40  const_cast<volScalarField&>
41  (
42  static_cast<const volScalarField&>
43  (
45  )
46  );
47 
49 
50  // Choose a master patch if one has not previously been chosen
51  if (masterPatchi_ == -1)
52  {
53  label masterPatchi = -1;
54 
55  forAll(bf, patchi)
56  {
57  if (isA<wallCellWallFunctionFvPatchScalarField>(bf[patchi]))
58  {
60  refCast<wallCellWallFunctionFvPatchScalarField>(bf[patchi]);
61 
62  if (masterPatchi == -1)
63  {
64  masterPatchi = patchi;
65  }
66 
67  epsilonPf.masterPatchi_ = masterPatchi;
68  }
69  }
70  }
71 
72  // Additional initialisation is only necessary on the master patch
73  if (patch().index() != masterPatchi_)
74  {
75  return;
76  }
77 
78  // Initialise at the start, and re-initialise following any mesh changes
79  if (wallCellsPtr_.valid() && !vf.mesh().changing())
80  {
81  return;
82  }
83 
84  // Construct flag lists for whether a cell or a (poly) boundary face is
85  // associated with a wall function
86  boolList cellIsWall(vf.mesh().nCells(), false);
87  boolList bFaceIsWall
88  (
89  vf.mesh().nFaces() - vf.mesh().nInternalFaces(),
90  false
91  );
92  forAll(bf, patchi)
93  {
94  if (isA<wallCellWallFunctionFvPatchScalarField>(bf[patchi]))
95  {
97  (
98  cellIsWall,
99  vf.mesh().boundary()[patchi].faceCells()
100  ) = true;
102  (
103  bFaceIsWall,
104  vf.mesh().polyFacesBf()[patchi] - vf.mesh().nInternalFaces()
105  ) = true;
106  }
107  }
108 
109  // Identify the wall cells
110  wallCellsPtr_.reset(new labelList(findIndices(cellIsWall, true)));
111 
112  // Construct the poly wall areas for each cell
113  scalarField cellMagWallArea(vf.mesh().nCells(), 0);
114  forAll(bFaceIsWall, bFacei)
115  {
116  if (bFaceIsWall[bFacei])
117  {
118  const label facei = bFacei + vf.mesh().nInternalFaces();
119  const label celli = vf.mesh().faceOwner()[facei];
120  cellMagWallArea[celli] += vf.mesh().magFaceAreas()[facei];
121  }
122  }
123 
124  // Construct the fv wall areas for each cell
125  scalarField cellMagSw(vf.mesh().nCells(), 0);
126  forAll(bf, patchi)
127  {
128  const fvPatch& fvp = vf.mesh().boundary()[patchi];
129 
130  if (isA<wallCellWallFunctionFvPatchScalarField>(bf[patchi]))
131  {
132  patchFieldAddToCellField(fvp, fvp.magSf(), cellMagSw);
133  }
134  }
135 
136  // Divide the areas to get the fraction of the cell values that is to be
137  // constrained. Modify this by the tolerance so that small fractions do not
138  // constrain anything.
139  tmp<scalarField> tWallCellFraction
140  (
141  scalarField(cellMagSw, wallCellsPtr_())
142  /scalarField(cellMagWallArea, wallCellsPtr_())
143  );
144  wallCellFractionPtr_.reset
145  (
146  max((tWallCellFraction - tol_)/(1 - tol_), scalar(0)).ptr()
147  );
148 }
149 
150 
152 (
153  const fvPatch& fvp,
154  const scalarField& pf,
155  scalarField& vf
156 )
157 {
158  const labelUList& patchFaceCells = fvp.faceCells();
159 
160  forAll(patchFaceCells, patchFacei)
161  {
162  vf[patchFaceCells[patchFacei]] += pf[patchFacei];
163  }
164 }
165 
166 
169 (
170  const PtrList<scalarField>& pfs
171 ) const
172 {
173  const volScalarField& vf =
174  static_cast<const volScalarField&>(internalField());
175 
176  const volScalarField::Boundary& bf = vf.boundaryField();
177 
178  scalarField cellMagSw(vf.mesh().nCells(), scalar(0));
179  scalarField cellMagSwValue(vf.mesh().nCells(), scalar(0));
180 
181  forAll(bf, patchi)
182  {
183  const fvPatch& fvp = vf.mesh().boundary()[patchi];
184 
185  if (isA<wallCellWallFunctionFvPatchScalarField>(bf[patchi]))
186  {
187  patchFieldAddToCellField
188  (
189  fvp,
190  fvp.magSf(),
191  cellMagSw
192  );
193  patchFieldAddToCellField
194  (
195  fvp,
196  fvp.magSf()*pfs[patchi],
197  cellMagSwValue
198  );
199  }
200  }
201 
202  return
203  scalarField(cellMagSwValue, wallCellsPtr_())
204  /max(scalarField(cellMagSw, wallCellsPtr_()), vSmall);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
209 
212 (
213  const fvPatch& p,
215  const dictionary& dict
216 )
217 :
218  fixedValueFvPatchField<scalar>(p, iF, dict, false),
219  masterPatchi_(-1),
220  wallCellsPtr_(nullptr),
221  wallCellFractionPtr_(nullptr)
222 {
223  // Apply a zero-gradient condition on start-up
225 }
226 
227 
230 (
232  const fvPatch& p,
234  const fieldMapper& mapper
235 )
236 :
237  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper, false),
238  masterPatchi_(-1),
239  wallCellsPtr_(nullptr),
240  wallCellFractionPtr_(nullptr)
241 {
242  // Apply a zero-gradient condition to any unmapped faces
243  mapper(*this, ptf, [&](){ return this->patchInternalField(); });
244 }
245 
246 
249 (
252 )
253 :
254  fixedValueFvPatchField<scalar>(ewfpsf, iF),
255  masterPatchi_(-1),
256  wallCellsPtr_(nullptr),
257  wallCellFractionPtr_(nullptr)
258 {}
259 
260 
261 
262 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
263 
265 (
266  const fvPatchScalarField& ptf,
267  const fieldMapper& mapper
268 )
269 {
270  // Apply a zero-gradient condition to any unmapped faces
271  mapper(*this, ptf, [&](){ return this->patchInternalField(); });
272  wallCellsPtr_.clear();
273  wallCellFractionPtr_.clear();
274 }
275 
276 
278 (
279  const fvPatchScalarField& ptf
280 )
281 {
283  wallCellsPtr_.clear();
284  wallCellFractionPtr_.clear();
285 }
286 
287 
288 // ************************************************************************* //
#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...
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A List with indirect addressing.
Definition: UIndirectList.H:60
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
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:362
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:415
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:170
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
Definition: fvPatchField.C:195
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:356
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:142
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
A class for managing temporary objects.
Definition: tmp.H:55
Base class for wall functions that modify cell values.
tmp< scalarField > patchFieldsToWallCellField(const PtrList< scalarField > &pfs) const
Average a set of patch fields into a wall cell field.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
static void patchFieldAddToCellField(const fvPatch &fvp, const scalarField &pf, scalarField &vf)
Sum values from a patch field into a cell field.
wallCellWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
void initMaster()
Initialise the master cell indices and fractions, and allocate the.
label patchi
const dimensionedScalar e
Elementary charge.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurrences of given element. Linear search.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict
volScalarField & p