Loading...
Searching...
No Matches
oversetAdjustPhi.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) 2017-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "oversetAdjustPhi.H"
29#include "volFields.H"
32#include "syncTools.H"
33#include "fv.H"
34
35// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
36
38(
40 const volVectorField& U,
41 const label zoneId
42)
43{
44 const fvMesh& mesh = U.mesh();
46 const labelUList& cellTypes = overlap.cellTypes();
47 const labelUList& zoneID = overlap.zoneID();
48
49 // Pass1: accumulate all fluxes, calculate correction factor
50
51 scalar massIn = 0;
52 scalar massOut = 0;
53
54 surfaceScalarField::Boundary& bphi = phi.boundaryFieldRef();
55
56 // Check all faces on the outside of interpolated cells
57 const labelUList& own = mesh.owner();
58 const labelUList& nei = mesh.neighbour();
59 forAll(own, facei)
60 {
61 const label zonei = zoneID[own[facei]];
62 const label ownType = cellTypes[own[facei]];
63 const label neiType = cellTypes[nei[facei]];
64
65 const bool ownCalc =
67 && (neiType == cellCellStencil::INTERPOLATED);
68
69 const bool neiCalc =
71 && (neiType == cellCellStencil::CALCULATED);
72
73 const bool ownOrCalc = (ownCalc || neiCalc);
74
75 if
76 (
77 (ownOrCalc && (zonei == zoneId))
78 || (ownOrCalc && (zoneId == -1))
79 )
80 {
81 // Calculate flux w.r.t. calculated cell
82 scalar flux = phi[facei];
83
84 if (ownCalc)
85 {
86 flux = -flux;
87 }
88
89 if (flux < 0.0)
90 {
91 massIn -= flux;
92 }
93 else
94 {
95 massOut += flux;
96 }
97 }
98 }
99
100
101 // Check all coupled faces on the outside of interpolated cells
102 labelList neiCellTypes;
104
105 forAll(bphi, patchi)
106 {
107 const fvPatchVectorField& Up = U.boundaryField()[patchi];
108 const fvsPatchScalarField& phip = bphi[patchi];
109 const labelUList& fc = Up.patch().faceCells();
110
111 const label start = Up.patch().start();
112
113 forAll(fc, i)
114 {
115 const label facei = start + i;
116 const label celli = fc[i];
117 const label ownType = cellTypes[celli];
118 const label neiType = neiCellTypes[facei - mesh.nInternalFaces()];
119
120 const label zonei = zoneID[celli];
121
122 const bool ownCalc =
123 (ownType == cellCellStencil::CALCULATED)
124 && (neiType == cellCellStencil::INTERPOLATED);
125
126 if (ownCalc && (zonei == zoneId))
127 {
128 // Calculate flux w.r.t. calculated cell
129 scalar flux = phip[i];
130
131 if (flux < 0.0)
132 {
133 massIn -= flux;
134 }
135 else
136 {
137 massOut += flux;
138 }
139 }
140 }
141 }
142 reduce(massIn, sumOp<scalar>());
143 reduce(massOut, sumOp<scalar>());
144
145 const scalar massCorr = massIn/(massOut + SMALL);
146
147 if (fv::debug)
148 {
149 Info<< "Zone : " << zoneId << nl
150 << "mass outflow : " << massOut << nl
151 << "mass inflow : " << massIn << nl
152 << "correction factor : " << massCorr << endl;
153 }
154
155
156 // Pass2: adjust fluxes
157 forAll(own, facei)
158 {
159 const label zonei = zoneID[own[facei]]; // own and nei in same zone
160
161 const label ownType = cellTypes[own[facei]];
162 const label neiType = cellTypes[nei[facei]];
163
164 const bool ownCalc =
165 (ownType == cellCellStencil::CALCULATED)
166 && (neiType == cellCellStencil::INTERPOLATED);
167
168 const bool neiCalc =
170 && (neiType == cellCellStencil::CALCULATED);
171
172 const bool ownOrCalc = (ownCalc || neiCalc);
173
174 if
175 (
176 (ownOrCalc && (zonei == zoneId)) || (ownOrCalc && (zoneId == -1))
177 )
178 {
179 scalar flux = phi[facei];
180
181 if (ownCalc)
182 {
183 flux = -flux;
184 }
185
186 if (flux < 0.0)
187 {
188 phi[facei] /= Foam::sqrt(massCorr);
189 }
190 else
191 {
192 phi[facei] *= Foam::sqrt(massCorr);
193 }
194 }
195 }
196
197 forAll(bphi, patchi)
198 {
199 const fvPatchVectorField& Up = U.boundaryField()[patchi];
200 fvsPatchScalarField& phip = bphi[patchi];
201 const labelUList& fc = Up.patch().faceCells();
202
203 const label start = Up.patch().start();
204
205 forAll(fc, i)
206 {
207 const label facei = start + i;
208 const label celli = fc[i];
209 const label zonei = zoneID[celli]; // note:own and nei in same zone
210 const label ownType = cellTypes[celli];
211 const label neiType = neiCellTypes[facei - mesh.nInternalFaces()];
212
213 const bool ownCalc =
214 (ownType == cellCellStencil::CALCULATED)
215 && (neiType == cellCellStencil::INTERPOLATED);
216
217 const bool neiCalc =
219 && (neiType == cellCellStencil::CALCULATED);
220
221 const bool ownOrCalc = (ownCalc || neiCalc);
222
223 if (ownOrCalc && (zonei == zoneId))
224 {
225 // Calculate flux w.r.t. calculated cell
226 scalar flux = phip[i];
227
228 if (neiCalc)
229 {
230 flux = -flux;
231 }
232
233 if (flux < 0.0)
234 {
235 // phip[i] /= Foam::sqrt(massCorr[zonei]);
236 }
237 else
238 {
239 phip[i] *= massCorr;
240 }
241 }
242 }
243 }
244
245
246 // Check correction
247 #ifdef FULLDEBUG
248 if (fv::debug)
249 {
250 scalar massOutCheck = 0;
251 scalar massInCheck = 0;
252
253 forAll(own, facei)
254 {
255 const label zonei = zoneID[own[facei]];
256 const label ownType = cellTypes[own[facei]];
257 const label neiType = cellTypes[nei[facei]];
258
259 const bool ownCalc =
260 (ownType == cellCellStencil::CALCULATED)
261 && (neiType == cellCellStencil::INTERPOLATED);
262
263 const bool neiCalc =
265 && (neiType == cellCellStencil::CALCULATED);
266
267 const bool ownOrCalc = (ownCalc || neiCalc);
268
269 if
270 (
271 (ownOrCalc && (zonei == zoneId))||(ownOrCalc && (zoneId == -1))
272 )
273 {
274 scalar flux = phi[facei];
275
276 if (ownCalc)
277 {
278 flux = -flux;
279 }
280
281 if (flux < 0.0)
282 {
283 massInCheck -= flux;
284 }
285 else
286 {
287 massOutCheck += flux;
288 }
289 }
290 }
291
292 Info<< "mass in:" << massInCheck << " out:" << massOutCheck << nl;
293 }
294 #endif /* FULLDEBUG */
295
296 return true;
297}
298
299
300// ************************************************************************* //
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
static FOAM_NO_DANGLING_REFERENCE const cellCellStencilObject & New(const fvMesh &mesh, Args &&... args)
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvPatch & patch() const noexcept
Return the patch.
label start() const noexcept
The patch start within the polyMesh face list.
Definition fvPatch.H:226
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData, const bool parRun=UPstream::parRun())
Extract and swap to obtain neighbour cell values for all boundary faces.
U
Definition pEqn.H:72
dynamicFvMesh & mesh
const cellCellStencilObject & overlap
Definition correctPhi.H:57
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool oversetAdjustPhi(surfaceScalarField &phi, const volVectorField &U, const label zoneId=-1)
Adjust the balance of fluxes to obey continuity.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvPatchField< vector > fvPatchVectorField
fvsPatchField< scalar > fvsPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity.
const labelUList & cellTypes
Definition setCellMask.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.