waWallFunctionFvPatchScalarField.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) 2014-2022 PCOpt/NTUA
9 Copyright (C) 2014-2022 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "wallFvPatch.H"
32#include "fvPatchFieldMapper.H"
33#include "volFields.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42
43void waWallFunctionFvPatchScalarField::checkType()
44{
45 if (!isA<wallFvPatch>(patch()))
46 {
48 << "Invalid wall function specification" << nl
49 << " Patch type for patch " << patch().name()
50 << " must be wall" << nl
51 << " Current patch type is " << patch().type() << nl << endl
52 << abort(FatalError);
53 }
54}
55
56
57// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58
60(
61 const fvPatch& p,
63)
64:
65 fixedValueFvPatchField<scalar>(p, iF),
67{
68 checkType();
69}
70
71
73(
75 const fvPatch& p,
77 const fvPatchFieldMapper& mapper
78)
79:
80 fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
81 adjointScalarBoundaryCondition(p, iF, ptf.adjointSolverName_)
82{
83 checkType();
84}
85
86
88(
89 const fvPatch& p,
91 const dictionary& dict
92)
93:
94 fixedValueFvPatchField<scalar>(p, iF, dict),
95 adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
96{
97 checkType();
98}
99
100
102(
104)
105:
106 fixedValueFvPatchField<scalar>(ewfpsf),
108{
109 checkType();
110}
111
112
114(
117)
118:
119 fixedValueFvPatchField<scalar>(ewfpsf, iF),
121{
122 checkType();
123}
124
125
126// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127
129(
130 fvMatrix<scalar>& matrix
131)
132{
133 scalarField& Diag = matrix.diag();
134 scalarField& lower = matrix.lower();
135 scalarField& upper = matrix.upper();
136 FieldField<Field, scalar>& internalCoeffs = matrix.internalCoeffs();
137 FieldField<Field, scalar>& boundaryCoeffs = matrix.boundaryCoeffs();
138 const fvMesh& mesh = patch().boundaryMesh().mesh();
139 const labelList& faceCells = patch().faceCells();
140
141 // Add diag term from the omega expression next to the wall
142 for (const label celli : faceCells)
143 {
144 Diag[celli] = 1;
145 }
146
147 // We want something similar to setValues, but slightly modified.
148 // The solution of the boundary cell should understand contributions from
149 // the second cells off the wall but they should see the
150 // solution of the boundary cell as zero.
151 // Contributions from neighbouring cells with an omegaWallFunction boundary
152 // condition should also be zero
153
154 const cellList& cells = mesh.cells();
155 const labelUList& own = mesh.owner();
156
157 /*
158 const labelUList& nei = mesh.neighbour();
159 const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
160 (
161 IOobject::groupName
162 (
163 turbulenceModel::propertiesName,
164 internalField().group()
165 )
166 );
167
168 const tmp<volScalarField> tomega = turbModel.omega();
169 const volScalarField& omega = tomega();
170 typedef omegaWallFunctionFvPatchScalarField omegaWF;
171 */
172
173 forAll(faceCells, i)
174 {
175 const label celli = faceCells[i];
176 const cell& c = cells[celli];
177
178 forAll(c, j)
179 {
180 const label facei = c[j];
181
182 if (mesh.isInternalFace(facei))
183 {
184 // Neighbouring cells should get no contribution from
185 // ourselves in all cases
186 //label cellNei(-1);
187 if (celli == own[facei])
188 {
189 //cellNei = nei[facei];
190 lower[facei] = 0.0;
191 }
192 else
193 {
194 //cellNei = own[facei];
195 upper[facei] = 0.0;
196 }
197 // Additionally, if the neighbouring cell is also a boundary
198 // one with omegaWallFunction in one of its faces,
199 // contributions between the two cells should be canceled out
200 // as well.
201 // Already covered by the above
202 /*
203 bool neiHasOmegaWFface(false);
204 const cell& neiCell = cells[cellNei];
205 forAll(neiCell, fNei)
206 {
207 const label faceNei = neiCell[fNei];
208
209 const label patchNei =
210 mesh.boundaryMesh().whichPatch(faceNei);
211 if (patchNei != -1)
212 {
213 const fvPatchField& omegaNei =
214 omega.boundaryField()[patchNei];
215 if (isA<omegaWF>(omegaNei))
216 {
217 neiHasOmegaWFface = true;
218 break;
219 }
220 }
221 }
222 if (neiHasOmegaWFface)
223 {
224 if (celli == own[facei])
225 {
226 upper[facei] = 0.0;
227 }
228 else
229 {
230 lower[facei] = 0.0;
231 }
232 }
233 */
234 }
235 // Contributions from boundaries should have already been removed
236 // using value*Coeffs and boundary*Coeffs
237 // Just to be safe
238 else
239 {
240 const label patchi = mesh.boundaryMesh().whichPatch(facei);
241 if (internalCoeffs[patchi].size())
242 {
243 label patchFacei =
244 mesh.boundaryMesh()[patchi].whichFace(facei);
245 internalCoeffs[patchi][patchFacei] = Zero;
246 boundaryCoeffs[patchi][patchFacei] = Zero;
247 }
248 }
249 }
250 }
251
253}
254
255
257{
258 if (updated())
259 {
260 return;
261 }
262
265}
266
267
269(
270 const tmp<scalarField>&
271) const
272{
273 return tmp<Field<scalar>>::New(this->size(), Zero);
274}
275
276
278(
279 const tmp<scalarField>&
280) const
281{
282 return tmp<Field<scalar>>::New(this->size(), Zero);
283}
284
285
288{
289 return tmp<Field<scalar>>::New(this->size(), Zero);
290}
291
292
295{
296 return tmp<Field<scalar>>::New(this->size(), Zero);
297}
298
299
301{
303}
304
305
306// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307
309(
312);
313
314// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315
316} // End namespace Foam
317
318// ************************************************************************* //
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...
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
label size() const noexcept
The number of elements in the UList.
Definition: UListI.H:420
Base class for solution control classes.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual bool write()
Write the output fields.
const fvMesh & mesh() const noexcept
Return the mesh reference.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
const FieldField< Field, Type > & internalCoeffs() const noexcept
Definition: fvMatrix.H:470
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
Definition: fvMatrix.H:484
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:417
A FieldMapper for finite-volume patch fields.
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:358
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:392
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual const word & name() const
Return name.
Definition: fvPatch.H:173
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:209
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:113
scalarField & upper()
Definition: lduMatrix.C:203
scalarField & lower()
Definition: lduMatrix.C:174
scalarField & diag()
Definition: lduMatrix.C:192
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
const cellList & cells() const
A class for managing temporary objects.
Definition: tmp.H:65
virtual tmp< Field< scalar > > gradientBoundaryCoeffs() const
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate waEqn at the cells next to the wall.
virtual tmp< Field< scalar > > valueBoundaryCoeffs(const tmp< scalarField > &) const
virtual tmp< Field< scalar > > valueInternalCoeffs(const tmp< scalarField > &) const
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< Field< scalar > > gradientInternalCoeffs() const
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
const cellShapeList & cells
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333