37template<
typename Type>
38Type Foam::isoAdvection::faceValue
40 const GeometricField<Type, fvsPatchField, surfaceMesh>&
f,
44 if (mesh_.isInternalFace(facei))
46 return f.primitiveField()[facei];
50 const polyBoundaryMesh&
pbm = mesh_.boundaryMesh();
53 const label patchi =
pbm.patchID(facei);
55 if (patchi < 0 || patchi >=
pbm.size())
58 <<
"Cannot find patch for face " << facei
63 const polyPatch&
pp =
pbm[patchi];
66 return pTraits<Type>::zero;
69 const label patchFacei =
pp.whichFace(facei);
70 return f.boundaryField()[patchi][patchFacei];
75template<
typename Type>
76void Foam::isoAdvection::setFaceValue
83 if (mesh_.isInternalFace(facei))
85 f.primitiveFieldRef()[facei] = value;
94 if (patchi < 0 || patchi >=
pbm.
size())
97 <<
"Cannot find patch for face " << facei
108 const label patchFacei =
pp.whichFace(facei);
109 f.boundaryFieldRef()[patchi][patchFacei] = value;
114template<
class SpType,
class SuType>
115void Foam::isoAdvection::limitFluxes
124 auto alpha1Limits =
gMinMax(alpha1_);
126 const scalar aTol = 100*SMALL;
127 scalar maxAlphaMinus1 = alpha1Limits.max() - 1;
128 scalar minAlpha = alpha1Limits.min();
129 const label nOvershoots = 20;
131 const labelList& owner = mesh_.faceOwner();
132 const labelList& neighbour = mesh_.faceNeighbour();
135 <<
"isoAdvection: Before conservative bounding: min(alpha) = "
136 << minAlpha <<
", max(alpha) = 1 + " << maxAlphaMinus1 <<
endl;
140 bitSet needBounding(mesh_.nCells(),
false);
141 needBounding.set(surfCells_);
143 extendMarkedCells(needBounding);
146 for (label
n = 0;
n < nAlphaBounds_;
n++)
148 if (maxAlphaMinus1 > aTol || minAlpha < -aTol)
153 dVfcorrectionValues =
Zero;
154 boundFlux(needBounding, dVfcorrectionValues, correctedFaces,
Sp,
Su);
156 correctedFaces.append
158 syncProcPatches(dVfcorrectionValues, phi_,
true)
162 for (
const label facei : correctedFaces)
164 if (alreadyUpdated.insert(facei))
166 checkIfOnProcPatch(facei);
167 const label own = owner[facei];
168 scalar Vown = mesh_.V()[own];
169 if (porosityEnabled_)
171 Vown *= porosityPtr_->primitiveField()[own];
173 alpha1_[own] -= faceValue(dVfcorrectionValues, facei)/Vown;
175 if (mesh_.isInternalFace(facei))
177 const label nei = neighbour[facei];
178 scalar Vnei = mesh_.V()[nei];
179 if (porosityEnabled_)
181 Vnei *= porosityPtr_->primitiveField()[nei];
184 faceValue(dVfcorrectionValues, facei)/Vnei;
188 const scalar corrVf =
189 faceValue(dVf_, facei)
190 + faceValue(dVfcorrectionValues, facei);
192 setFaceValue(dVf_, facei, corrVf);
196 if (!mesh_.isInternalFace(facei))
199 mesh_.boundaryMesh();
203 const label patchFacei =
pp.whichFace(facei);
207 const label neiPatchID(cpp->neighbPolyPatchID());
209 dVf_.boundaryFieldRef();
210 dVfb[neiPatchID][patchFacei] =
211 -dVfb[patchi][patchFacei];
213 const label globalFacei = np.start() + patchFacei;
214 const label neiOwn(owner[globalFacei]);
215 scalar VneiOwn = mesh_.V()[neiOwn];
216 if (porosityEnabled_)
219 porosityPtr_->primitiveField()[neiOwn];
222 faceValue(dVfcorrectionValues, facei)/VneiOwn;
230 syncProcPatches(dVf_, phi_);
237 alpha1Limits =
gMinMax(alpha1_);
239 maxAlphaMinus1 = alpha1Limits.max() - 1;
240 minAlpha = alpha1Limits.min();
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));
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;
259 alpha1_.correctBoundaryConditions();
264template<
class SpType,
class SuType>
265void Foam::isoAdvection::boundFlux
267 const bitSet& nextToInterface,
276 const scalar rDeltaT = 1.0/mesh_.time().deltaTValue();
278 correctedFaces.
clear();
279 const scalar aTol = 100*SMALL;
282 const scalar dt = mesh_.time().deltaTValue();
291 for (label celli: nextToInterface)
293 if (alpha1_[celli] < -aTol || alpha1_[celli] > 1 + aTol)
295 scalar Vi = meshV[celli];
296 if (porosityEnabled_)
298 Vi *= porosityPtr_->primitiveField()[celli];
301 scalar alphaOvershoot =
302 pos0(alpha1_[celli] - 1)*(alpha1_[celli] - 1)
303 +
neg0(alpha1_[celli])*alpha1_[celli];
305 scalar fluidToPassOn = alphaOvershoot*Vi;
306 label nFacesToPassFluidThrough = 1;
308 bool firstLoop =
true;
312 for (label iter=0; iter<10; iter++)
314 if (
mag(alphaOvershoot) < aTol || nFacesToPassFluidThrough == 0)
320 <<
"\n\nBounding cell " << celli
321 <<
" with alpha overshooting " << alphaOvershoot
324 facesToPassFluidThrough.
clear();
329 setDownwindFaces(celli, downwindFaces);
332 nFacesToPassFluidThrough = 0;
334 for (
const label facei : downwindFaces)
336 const scalar phif = faceValue(phi_, facei);
339 faceValue(dVf_, facei)
340 + faceValue(dVfcorrectionValues, facei);
342 const scalar maxExtraFaceFluidTrans =
343 mag(
pos0(fluidToPassOn)*phif*dt - dVff);
352 <<
"downwindFace " << facei
353 <<
" has maxExtraFaceFluidTrans = "
354 << maxExtraFaceFluidTrans <<
endl;
356 if (maxExtraFaceFluidTrans/Vi > aTol)
362 facesToPassFluidThrough.
append(facei);
364 dVfmax.
append(maxExtraFaceFluidTrans);
365 dVftot +=
mag(phif*dt);
370 <<
"\nfacesToPassFluidThrough: "
371 << facesToPassFluidThrough <<
", dVftot = "
372 << dVftot <<
" m3 corresponding to dalpha = "
373 << dVftot/Vi <<
endl;
375 forAll(facesToPassFluidThrough, fi)
377 const label facei = facesToPassFluidThrough[fi];
378 scalar fluidToPassThroughFace =
379 mag(fluidToPassOn)*
mag(
phi[fi]*dt)/dVftot;
381 nFacesToPassFluidThrough +=
382 pos0(dVfmax[fi] - fluidToPassThroughFace);
384 fluidToPassThroughFace =
385 min(fluidToPassThroughFace, dVfmax[fi]);
387 scalar dVff = faceValue(dVfcorrectionValues, facei);
390 sign(
phi[fi])*fluidToPassThroughFace*
sign(fluidToPassOn);
392 setFaceValue(dVfcorrectionValues, facei, dVff);
396 checkIfOnProcPatch(facei);
397 correctedFaces.
append(facei);
405 alphaOld[celli]*rDeltaT +
Su[celli]
406 - netFlux(dVf_, celli)/Vi*rDeltaT
407 - netFlux(dVfcorrectionValues, celli)/Vi*rDeltaT
410 (rDeltaT -
Sp[celli]);
413 pos0(alpha1New - 1)*(alpha1New - 1)
414 +
neg0(alpha1New)*alpha1New;
416 fluidToPassOn = alphaOvershoot*Vi;
419 <<
"\nNew alpha for cell " << celli <<
": "
420 << alpha1New <<
endl;
429template<
class SpType,
class SuType>
435 if (mesh_.topoChanging())
437 setProcessorPatches();
440 scalar advectionStartTime = mesh_.time().elapsedCpuTime();
442 const scalar rDeltaT = 1.0/mesh_.time().deltaTValue();
445 surf_->reconstruct();
447 if (timeIndex_ < mesh_.time().timeIndex())
449 timeIndex_= mesh_.time().timeIndex();
450 surf_->normal().oldTime() = surf_->normal();
451 surf_->centre().oldTime() = surf_->centre();
459 timeIntegratedFlux();
464 alpha1In_ *= (mesh_.Vsc0()/mesh_.Vsc());
468 if (porosityEnabled_)
471 alpha1_.primitiveFieldRef() =
473 alpha1_.oldTime().primitiveField()*rDeltaT +
Su.field()
475 /porosityPtr_->primitiveField()
476 )/(rDeltaT -
Sp.field());
480 alpha1_.primitiveFieldRef() =
482 alpha1_.oldTime().primitiveField()*rDeltaT
485 )/(rDeltaT -
Sp.field());
488 alpha1_.correctBoundaryConditions();
496 Info<<
"isoAdvection: After conservative bounding:"
497 <<
" min(alpha) = " <<
limits.min()
498 <<
", max(alpha) = 1 + " << (
limits.max()-1) <<
endl;
503 applyBruteForceBounding();
508 advectionTime_ += (mesh_.time().elapsedCpuTime() - advectionStartTime);
510 <<
"isoAdvection: time consumption = "
511 << label(100*advectionTime_/(mesh_.time().elapsedCpuTime() + SMALL))
514 alphaPhi_ = dVf_/mesh_.time().deltaT();
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.
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.
label size() const noexcept
The number of entries in the list.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
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.
Upwind differencing scheme class.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
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.
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
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.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
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.
errorManip< error > abort(error &err)
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).
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).
#define addProfilingInFunction(Name)
Define profiling trigger with specified name and description corresponding to the compiler-defined fu...
#define forAll(list, i)
Loop across all elements in list.