Loading...
Searching...
No Matches
isoAdvectionTemplates.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) 2016-2017 DHI
9 Modified code Copyright (C) 2016-2017 OpenCFD Ltd.
10 Modified code Copyright (C) 2019-2020 DLR
11 Modified code Copyright (C) 2021 Johan Roenby
12-------------------------------------------------------------------------------
13License
14 This file is part of OpenFOAM.
15
16 OpenFOAM is free software: you can redistribute it and/or modify it
17 under the terms of the GNU General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24 for more details.
25
26 You should have received a copy of the GNU General Public License
27 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
31#include "isoAdvection.H"
32#include "fvcSurfaceIntegrate.H"
33#include "upwind.H"
34
35// ************************************************************************* //
36
37template<typename Type>
38Type Foam::isoAdvection::faceValue
39(
40 const GeometricField<Type, fvsPatchField, surfaceMesh>& f,
41 const label facei
42) const
43{
44 if (mesh_.isInternalFace(facei))
45 {
46 return f.primitiveField()[facei];
47 }
48 else
49 {
50 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
51
52 // Boundary face. Find out which face of which patch
53 const label patchi = pbm.patchID(facei);
54
55 if (patchi < 0 || patchi >= pbm.size())
56 {
58 << "Cannot find patch for face " << facei
59 << abort(FatalError);
60 }
61
62 // Handle empty patches
63 const polyPatch& pp = pbm[patchi];
64 if (isA<emptyPolyPatch>(pp) || pp.empty())
65 {
66 return pTraits<Type>::zero;
67 }
68
69 const label patchFacei = pp.whichFace(facei);
70 return f.boundaryField()[patchi][patchFacei];
71 }
72}
73
74
75template<typename Type>
76void Foam::isoAdvection::setFaceValue
77(
79 const label facei,
80 const Type& value
81) const
82{
83 if (mesh_.isInternalFace(facei))
84 {
85 f.primitiveFieldRef()[facei] = value;
86 }
87 else
88 {
89 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
90
91 // Boundary face. Find out which face of which patch
92 const label patchi = pbm.patchID(facei);
93
94 if (patchi < 0 || patchi >= pbm.size())
95 {
97 << "Cannot find patch for face " << facei
98 << abort(FatalError);
99 }
100
101 // Handle empty patches
102 const polyPatch& pp = pbm[patchi];
103 if (isA<emptyPolyPatch>(pp) || pp.empty())
104 {
105 return;
106 }
107
108 const label patchFacei = pp.whichFace(facei);
109 f.boundaryFieldRef()[patchi][patchFacei] = value;
110 }
111}
112
113
114template<class SpType, class SuType>
115void Foam::isoAdvection::limitFluxes
116(
117 const SpType& Sp,
118 const SuType& Su
119)
120{
121 addProfilingInFunction(geometricVoF);
123
124 auto alpha1Limits = gMinMax(alpha1_);
125
126 const scalar aTol = 100*SMALL; // Note: tolerances
127 scalar maxAlphaMinus1 = alpha1Limits.max() - 1; // max(alphaNew - 1);
128 scalar minAlpha = alpha1Limits.min(); // min(alphaNew);
129 const label nOvershoots = 20; // sum(pos0(alphaNew - 1 - aTol));
130
131 const labelList& owner = mesh_.faceOwner();
132 const labelList& neighbour = mesh_.faceNeighbour();
133
135 << "isoAdvection: Before conservative bounding: min(alpha) = "
136 << minAlpha << ", max(alpha) = 1 + " << maxAlphaMinus1 << endl;
137
138 surfaceScalarField dVfcorrectionValues("dVfcorrectionValues", dVf_*0.0);
139
140 bitSet needBounding(mesh_.nCells(),false);
141 needBounding.set(surfCells_);
142
143 extendMarkedCells(needBounding);
144
145 // Loop number of bounding steps
146 for (label n = 0; n < nAlphaBounds_; n++)
147 {
148 if (maxAlphaMinus1 > aTol || minAlpha < -aTol) // Note: tolerances
149 {
150 DebugInfo << "boundAlpha... " << endl;
151
152 DynamicList<label> correctedFaces(3*nOvershoots);
153 dVfcorrectionValues = Zero;
154 boundFlux(needBounding, dVfcorrectionValues, correctedFaces,Sp,Su);
155
156 correctedFaces.append
157 (
158 syncProcPatches(dVfcorrectionValues, phi_,true)
159 );
160
161 labelHashSet alreadyUpdated;
162 for (const label facei : correctedFaces)
163 {
164 if (alreadyUpdated.insert(facei))
165 {
166 checkIfOnProcPatch(facei);
167 const label own = owner[facei];
168 scalar Vown = mesh_.V()[own];
169 if (porosityEnabled_)
170 {
171 Vown *= porosityPtr_->primitiveField()[own];
172 }
173 alpha1_[own] -= faceValue(dVfcorrectionValues, facei)/Vown;
174
175 if (mesh_.isInternalFace(facei))
176 {
177 const label nei = neighbour[facei];
178 scalar Vnei = mesh_.V()[nei];
179 if (porosityEnabled_)
180 {
181 Vnei *= porosityPtr_->primitiveField()[nei];
182 }
183 alpha1_[nei] +=
184 faceValue(dVfcorrectionValues, facei)/Vnei;
185 }
186
187 // Change to treat boundaries consistently
188 const scalar corrVf =
189 faceValue(dVf_, facei)
190 + faceValue(dVfcorrectionValues, facei);
191
192 setFaceValue(dVf_, facei, corrVf);
193
194 // If facei is on a cyclic patch correct dVf of facej on
195 // neighbour patch and alpha value of facej's owner cell.
196 if (!mesh_.isInternalFace(facei))
197 {
199 mesh_.boundaryMesh();
200 const label patchi = boundaryMesh.patchID(facei);
201 const polyPatch& pp = boundaryMesh[patchi];
203 const label patchFacei = pp.whichFace(facei);
204
205 if (cpp)
206 {
207 const label neiPatchID(cpp->neighbPolyPatchID());
209 dVf_.boundaryFieldRef();
210 dVfb[neiPatchID][patchFacei] =
211 -dVfb[patchi][patchFacei];
212 const polyPatch& np = boundaryMesh[neiPatchID];
213 const label globalFacei = np.start() + patchFacei;
214 const label neiOwn(owner[globalFacei]);
215 scalar VneiOwn = mesh_.V()[neiOwn];
216 if (porosityEnabled_)
217 {
218 VneiOwn *=
219 porosityPtr_->primitiveField()[neiOwn];
220 }
221 alpha1_[neiOwn] +=
222 faceValue(dVfcorrectionValues, facei)/VneiOwn;
223 }
224 }
225
226 }
227
228 }
229
230 syncProcPatches(dVf_, phi_);
231 }
232 else
233 {
234 break;
235 }
236
237 alpha1Limits = gMinMax(alpha1_);
238
239 maxAlphaMinus1 = alpha1Limits.max() - 1; // max(alphaNew - 1);
240 minAlpha = alpha1Limits.min(); // min(alphaNew);
241
242 if (debug)
243 {
244 // Check if still unbounded
245 //scalarField alphaNew(alpha1In_ - fvc::surfaceIntegrate(dVf_)());
246 label maxAlphaMinus1 = max(alpha1_.primitiveField() - 1);
247 scalar minAlpha = min(alpha1_.primitiveField());
248 label nUndershoots = sum(neg0(alpha1_.primitiveField() + aTol));
249 label nOvershoots = sum(pos0(alpha1_.primitiveField() - 1 - aTol));
250
251 Info<< "After bounding number " << n + 1 << " of time "
252 << mesh_.time().value() << ":" << nl
253 << "nOvershoots = " << nOvershoots << " with max(alpha1_-1) = "
254 << maxAlphaMinus1 << " and nUndershoots = " << nUndershoots
255 << " with min(alpha1_) = " << minAlpha << endl;
256 }
257 }
259 alpha1_.correctBoundaryConditions();
260
261}
262
263
264template<class SpType, class SuType>
265void Foam::isoAdvection::boundFlux
266(
267 const bitSet& nextToInterface,
268 surfaceScalarField& dVfcorrectionValues,
269 DynamicList<label>& correctedFaces,
270 const SpType& Sp,
271 const SuType& Su
272)
273{
274 addProfilingInFunction(geometricVoF);
276 const scalar rDeltaT = 1.0/mesh_.time().deltaTValue();
277
278 correctedFaces.clear();
279 const scalar aTol = 100*SMALL; // Note: tolerances
280
281 const scalarField& meshV = mesh_.cellVolumes();
282 const scalar dt = mesh_.time().deltaTValue();
283
284 DynamicList<label> downwindFaces(10);
285 DynamicList<label> facesToPassFluidThrough(downwindFaces.size());
286 DynamicList<scalar> dVfmax(downwindFaces.size());
287 DynamicList<scalar> phi(downwindFaces.size());
288 const volScalarField& alphaOld = alpha1_.oldTime();
289
290 // Loop through alpha cell centred field
291 for (label celli: nextToInterface)
292 {
293 if (alpha1_[celli] < -aTol || alpha1_[celli] > 1 + aTol)
294 {
295 scalar Vi = meshV[celli];
296 if (porosityEnabled_)
297 {
298 Vi *= porosityPtr_->primitiveField()[celli];
299 }
300
301 scalar alphaOvershoot =
302 pos0(alpha1_[celli] - 1)*(alpha1_[celli] - 1)
303 + neg0(alpha1_[celli])*alpha1_[celli];
304
305 scalar fluidToPassOn = alphaOvershoot*Vi;
306 label nFacesToPassFluidThrough = 1;
307
308 bool firstLoop = true;
309
310 // First try to pass surplus fluid on to neighbour cells that are
311 // not filled and to which dVf < phi*dt
312 for (label iter=0; iter<10; iter++)
313 {
314 if (mag(alphaOvershoot) < aTol || nFacesToPassFluidThrough == 0)
315 {
316 break;
317 }
318
320 << "\n\nBounding cell " << celli
321 << " with alpha overshooting " << alphaOvershoot
322 << endl;
323
324 facesToPassFluidThrough.clear();
325 dVfmax.clear();
326 phi.clear();
327
328 // Find potential neighbour cells to pass surplus phase to
329 setDownwindFaces(celli, downwindFaces);
330
331 scalar dVftot = 0;
332 nFacesToPassFluidThrough = 0;
333
334 for (const label facei : downwindFaces)
335 {
336 const scalar phif = faceValue(phi_, facei);
337
338 const scalar dVff =
339 faceValue(dVf_, facei)
340 + faceValue(dVfcorrectionValues, facei);
341
342 const scalar maxExtraFaceFluidTrans =
343 mag(pos0(fluidToPassOn)*phif*dt - dVff);
344
345 // dVf has same sign as phi and so if phi>0 we have
346 // mag(phi_[facei]*dt) - mag(dVf[facei]) = phi_[facei]*dt
347 // - dVf[facei]
348 // If phi < 0 we have mag(phi_[facei]*dt) -
349 // mag(dVf[facei]) = -phi_[facei]*dt - (-dVf[facei]) > 0
350 // since mag(dVf) < phi*dt
352 << "downwindFace " << facei
353 << " has maxExtraFaceFluidTrans = "
354 << maxExtraFaceFluidTrans << endl;
355
356 if (maxExtraFaceFluidTrans/Vi > aTol)
357 {
358// if (maxExtraFaceFluidTrans/Vi > aTol &&
359// mag(dVfIn[facei])/Vi > aTol) //Last condition may be
360// important because without this we will flux through uncut
361// downwind faces
362 facesToPassFluidThrough.append(facei);
363 phi.append(phif);
364 dVfmax.append(maxExtraFaceFluidTrans);
365 dVftot += mag(phif*dt);
366 }
367 }
368
370 << "\nfacesToPassFluidThrough: "
371 << facesToPassFluidThrough << ", dVftot = "
372 << dVftot << " m3 corresponding to dalpha = "
373 << dVftot/Vi << endl;
374
375 forAll(facesToPassFluidThrough, fi)
376 {
377 const label facei = facesToPassFluidThrough[fi];
378 scalar fluidToPassThroughFace =
379 mag(fluidToPassOn)*mag(phi[fi]*dt)/dVftot;
380
381 nFacesToPassFluidThrough +=
382 pos0(dVfmax[fi] - fluidToPassThroughFace);
383
384 fluidToPassThroughFace =
385 min(fluidToPassThroughFace, dVfmax[fi]);
386
387 scalar dVff = faceValue(dVfcorrectionValues, facei);
388
389 dVff +=
390 sign(phi[fi])*fluidToPassThroughFace*sign(fluidToPassOn);
391
392 setFaceValue(dVfcorrectionValues, facei, dVff);
393
394 if (firstLoop)
395 {
396 checkIfOnProcPatch(facei);
397 correctedFaces.append(facei);
398 }
399 }
400
401 firstLoop = false;
402
403 scalar alpha1New =
404 (
405 alphaOld[celli]*rDeltaT + Su[celli]
406 - netFlux(dVf_, celli)/Vi*rDeltaT
407 - netFlux(dVfcorrectionValues, celli)/Vi*rDeltaT
408 )
409 /
410 (rDeltaT - Sp[celli]);
411
412 alphaOvershoot =
413 pos0(alpha1New - 1)*(alpha1New - 1)
414 + neg0(alpha1New)*alpha1New;
415
416 fluidToPassOn = alphaOvershoot*Vi;
417
419 << "\nNew alpha for cell " << celli << ": "
420 << alpha1New << endl;
421 }
422 }
424
425 DebugInfo << "correctedFaces = " << correctedFaces << endl;
426}
427
428
429template<class SpType, class SuType>
430void Foam::isoAdvection::advect(const SpType& Sp, const SuType& Su)
431{
432 addProfilingInFunction(geometricVoF);
434
435 if (mesh_.topoChanging())
436 {
437 setProcessorPatches();
438 }
439
440 scalar advectionStartTime = mesh_.time().elapsedCpuTime();
441
442 const scalar rDeltaT = 1.0/mesh_.time().deltaTValue();
443
444 // reconstruct the interface
445 surf_->reconstruct();
446
447 if (timeIndex_ < mesh_.time().timeIndex())
448 {
449 timeIndex_= mesh_.time().timeIndex();
450 surf_->normal().oldTime() = surf_->normal();
451 surf_->centre().oldTime() = surf_->centre();
452 }
453
454 // Initialising dVf with upwind values
455 // i.e. phi[facei]*alpha1[upwindCell[facei]]*dt
456 dVf_ = upwind<scalar>(mesh_, phi_).flux(alpha1_)*mesh_.time().deltaT();
457
458 // Do the isoAdvection on surface cells
459 timeIntegratedFlux();
460
461 // Adjust alpha for mesh motion
462 if (mesh_.moving())
463 {
464 alpha1In_ *= (mesh_.Vsc0()/mesh_.Vsc());
465 }
466
467 // Advect the free surface
468 if (porosityEnabled_)
469 {
470 // Should Su and Sp also be divided by porosity?
471 alpha1_.primitiveFieldRef() =
472 (
473 alpha1_.oldTime().primitiveField()*rDeltaT + Su.field()
474 - fvc::surfaceIntegrate(dVf_)().primitiveField()*rDeltaT
475 /porosityPtr_->primitiveField()
476 )/(rDeltaT - Sp.field());
477 }
478 else
479 {
480 alpha1_.primitiveFieldRef() =
481 (
482 alpha1_.oldTime().primitiveField()*rDeltaT
483 + Su.field()
484 - fvc::surfaceIntegrate(dVf_)().primitiveField()*rDeltaT
485 )/(rDeltaT - Sp.field());
486 }
487
488 alpha1_.correctBoundaryConditions();
489
490 // Adjust dVf for unbounded cells
491 limitFluxes(Sp, Su);
492
493 {
494 auto limits = gMinMax(alpha1In_);
495
496 Info<< "isoAdvection: After conservative bounding:"
497 << " min(alpha) = " << limits.min()
498 << ", max(alpha) = 1 + " << (limits.max()-1) << endl;
499 }
500
501 // Apply non-conservative bounding mechanisms (clipping and snapping)
502 // Note: We should be able to write out alpha before this is done!
503 applyBruteForceBounding();
504
505 // Write surface cell set and bound cell set if required by user
506 writeSurfaceCells();
507
508 advectionTime_ += (mesh_.time().elapsedCpuTime() - advectionStartTime);
510 << "isoAdvection: time consumption = "
511 << label(100*advectionTime_/(mesh_.time().elapsedCpuTime() + SMALL))
512 << "%" << endl;
513
514 alphaPhi_ = dVf_/mesh_.time().deltaT();
515}
516
517
518// ************************************************************************* //
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
Generic GeometricField class.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
bool empty() const noexcept
True if the list is empty (ie, size() is zero).
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Cyclic plane patch.
void writeSurfaceCells() const
Return cellSet of surface cells.
void advect(const SpType &Sp, const SuType &Su)
Advect the free surface. Updates alpha field, taking into account.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const labelList & patchID() const
Per boundary face label the patch index.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Upwind differencing scheme class.
Definition upwind.H:55
auto limits
Definition setRDeltaT.H:186
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
zeroField Su
Definition alphaSuSp.H:1
zeroField Sp
Definition alphaSuSp.H:2
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar neg0(const dimensionedScalar &ds)
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addProfilingInFunction(Name)
Define profiling trigger with specified name and description corresponding to the compiler-defined fu...
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299