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 -------------------------------------------------------------------------------
13 License
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 
37 template<typename Type>
38 Type 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 
75 template<typename Type>
76 void Foam::isoAdvection::setFaceValue
77 (
78  GeometricField<Type, fvsPatchField, surfaceMesh>& f,
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 
114 template<class SpType, class SuType>
115 void Foam::isoAdvection::limitFluxes
116 (
117  const SpType& Sp,
118  const SuType& Su
119 )
120 {
121  addProfilingInFunction(geometricVoF);
123 
124  const scalar aTol = 100*SMALL; // Note: tolerances
125  scalar maxAlphaMinus1 = gMax(alpha1_) - 1; // max(alphaNew - 1);
126  scalar minAlpha = gMin(alpha1_); // min(alphaNew);
127  const label nOvershoots = 20; // sum(pos0(alphaNew - 1 - aTol));
128 
129  const labelList& owner = mesh_.faceOwner();
130  const labelList& neighbour = mesh_.faceNeighbour();
131 
132  DebugInfo
133  << "isoAdvection: Before conservative bounding: min(alpha) = "
134  << minAlpha << ", max(alpha) = 1 + " << maxAlphaMinus1 << endl;
135 
136  surfaceScalarField dVfcorrectionValues("dVfcorrectionValues", dVf_*0.0);
137 
138  bitSet needBounding(mesh_.nCells(),false);
139  needBounding.set(surfCells_);
140 
141  extendMarkedCells(needBounding);
142 
143  // Loop number of bounding steps
144  for (label n = 0; n < nAlphaBounds_; n++)
145  {
146  if (maxAlphaMinus1 > aTol || minAlpha < -aTol) // Note: tolerances
147  {
148  DebugInfo << "boundAlpha... " << endl;
149 
150  DynamicList<label> correctedFaces(3*nOvershoots);
151  dVfcorrectionValues = dimensionedScalar("0",dimVolume,0.0);
152  boundFlux(needBounding, dVfcorrectionValues, correctedFaces,Sp,Su);
153 
154  correctedFaces.append
155  (
156  syncProcPatches(dVfcorrectionValues, phi_,true)
157  );
158 
159  labelHashSet alreadyUpdated;
160  for (const label facei : correctedFaces)
161  {
162  if (alreadyUpdated.insert(facei))
163  {
164  checkIfOnProcPatch(facei);
165  const label own = owner[facei];
166  scalar Vown = mesh_.V()[own];
167  if (porosityEnabled_)
168  {
169  Vown *= porosityPtr_->primitiveField()[own];
170  }
171  alpha1_[own] -= faceValue(dVfcorrectionValues, facei)/Vown;
172 
173  if (mesh_.isInternalFace(facei))
174  {
175  const label nei = neighbour[facei];
176  scalar Vnei = mesh_.V()[nei];
177  if (porosityEnabled_)
178  {
179  Vnei *= porosityPtr_->primitiveField()[nei];
180  }
181  alpha1_[nei] +=
182  faceValue(dVfcorrectionValues, facei)/Vnei;
183  }
184 
185  // Change to treat boundaries consistently
186  const scalar corrVf =
187  faceValue(dVf_, facei)
188  + faceValue(dVfcorrectionValues, facei);
189 
190  setFaceValue(dVf_, facei, corrVf);
191  }
192  }
193  syncProcPatches(dVf_, phi_);
194  }
195  else
196  {
197  break;
198  }
199 
200  maxAlphaMinus1 = gMax(alpha1_) - 1; // max(alphaNew - 1);
201  minAlpha = gMin(alpha1_); // min(alphaNew);
202 
203  if (debug)
204  {
205  // Check if still unbounded
206  //scalarField alphaNew(alpha1In_ - fvc::surfaceIntegrate(dVf_)());
207  label maxAlphaMinus1 = max(alpha1_.primitiveField() - 1);
208  scalar minAlpha = min(alpha1_.primitiveField());
209  label nUndershoots = sum(neg0(alpha1_.primitiveField() + aTol));
210  label nOvershoots = sum(pos0(alpha1_.primitiveField() - 1 - aTol));
211 
212  Info<< "After bounding number " << n + 1 << " of time "
213  << mesh_.time().value() << ":" << nl
214  << "nOvershoots = " << nOvershoots << " with max(alpha1_-1) = "
215  << maxAlphaMinus1 << " and nUndershoots = " << nUndershoots
216  << " with min(alpha1_) = " << minAlpha << endl;
217  }
218  }
219 
220  alpha1_.correctBoundaryConditions();
221 
222 }
223 
224 
225 template<class SpType, class SuType>
226 void Foam::isoAdvection::boundFlux
227 (
228  const bitSet& nextToInterface,
229  surfaceScalarField& dVfcorrectionValues,
230  DynamicList<label>& correctedFaces,
231  const SpType& Sp,
232  const SuType& Su
233 )
234 {
235  addProfilingInFunction(geometricVoF);
237  const scalar rDeltaT = 1.0/mesh_.time().deltaTValue();
238 
239  correctedFaces.clear();
240  const scalar aTol = 100*SMALL; // Note: tolerances
241 
242  const scalarField& meshV = mesh_.cellVolumes();
243  const scalar dt = mesh_.time().deltaTValue();
244 
245  DynamicList<label> downwindFaces(10);
246  DynamicList<label> facesToPassFluidThrough(downwindFaces.size());
247  DynamicList<scalar> dVfmax(downwindFaces.size());
248  DynamicList<scalar> phi(downwindFaces.size());
249  const volScalarField& alphaOld = alpha1_.oldTime();
250 
251  // Loop through alpha cell centred field
252  for (label celli: nextToInterface)
253  {
254  if (alpha1_[celli] < -aTol || alpha1_[celli] > 1 + aTol)
255  {
256  scalar Vi = meshV[celli];
257  if (porosityEnabled_)
258  {
259  Vi *= porosityPtr_->primitiveField()[celli];
260  }
261 
262  scalar alphaOvershoot =
263  pos0(alpha1_[celli] - 1)*(alpha1_[celli] - 1)
264  + neg0(alpha1_[celli])*alpha1_[celli];
265 
266  scalar fluidToPassOn = alphaOvershoot*Vi;
267  label nFacesToPassFluidThrough = 1;
268 
269  bool firstLoop = true;
270 
271  // First try to pass surplus fluid on to neighbour cells that are
272  // not filled and to which dVf < phi*dt
273  for (label iter=0; iter<10; iter++)
274  {
275  if (mag(alphaOvershoot) < aTol || nFacesToPassFluidThrough == 0)
276  {
277  break;
278  }
279 
280  DebugInfo
281  << "\n\nBounding cell " << celli
282  << " with alpha overshooting " << alphaOvershoot
283  << endl;
284 
285  facesToPassFluidThrough.clear();
286  dVfmax.clear();
287  phi.clear();
288 
289  // Find potential neighbour cells to pass surplus phase to
290  setDownwindFaces(celli, downwindFaces);
291 
292  scalar dVftot = 0;
293  nFacesToPassFluidThrough = 0;
294 
295  for (const label facei : downwindFaces)
296  {
297  const scalar phif = faceValue(phi_, facei);
298 
299  const scalar dVff =
300  faceValue(dVf_, facei)
301  + faceValue(dVfcorrectionValues, facei);
302 
303  const scalar maxExtraFaceFluidTrans =
304  mag(pos0(fluidToPassOn)*phif*dt - dVff);
305 
306  // dVf has same sign as phi and so if phi>0 we have
307  // mag(phi_[facei]*dt) - mag(dVf[facei]) = phi_[facei]*dt
308  // - dVf[facei]
309  // If phi < 0 we have mag(phi_[facei]*dt) -
310  // mag(dVf[facei]) = -phi_[facei]*dt - (-dVf[facei]) > 0
311  // since mag(dVf) < phi*dt
312  DebugInfo
313  << "downwindFace " << facei
314  << " has maxExtraFaceFluidTrans = "
315  << maxExtraFaceFluidTrans << endl;
316 
317  if (maxExtraFaceFluidTrans/Vi > aTol)
318  {
319 // if (maxExtraFaceFluidTrans/Vi > aTol &&
320 // mag(dVfIn[facei])/Vi > aTol) //Last condition may be
321 // important because without this we will flux through uncut
322 // downwind faces
323  facesToPassFluidThrough.append(facei);
324  phi.append(phif);
325  dVfmax.append(maxExtraFaceFluidTrans);
326  dVftot += mag(phif*dt);
327  }
328  }
329 
330  DebugInfo
331  << "\nfacesToPassFluidThrough: "
332  << facesToPassFluidThrough << ", dVftot = "
333  << dVftot << " m3 corresponding to dalpha = "
334  << dVftot/Vi << endl;
335 
336  forAll(facesToPassFluidThrough, fi)
337  {
338  const label facei = facesToPassFluidThrough[fi];
339  scalar fluidToPassThroughFace =
340  mag(fluidToPassOn)*mag(phi[fi]*dt)/dVftot;
341 
342  nFacesToPassFluidThrough +=
343  pos0(dVfmax[fi] - fluidToPassThroughFace);
344 
345  fluidToPassThroughFace =
346  min(fluidToPassThroughFace, dVfmax[fi]);
347 
348  scalar dVff = faceValue(dVfcorrectionValues, facei);
349 
350  dVff +=
351  sign(phi[fi])*fluidToPassThroughFace*sign(fluidToPassOn);
352 
353  setFaceValue(dVfcorrectionValues, facei, dVff);
354 
355  if (firstLoop)
356  {
357  checkIfOnProcPatch(facei);
358  correctedFaces.append(facei);
359  }
360  }
361 
362  firstLoop = false;
363 
364  scalar alpha1New =
365  (
366  alphaOld[celli]*rDeltaT + Su[celli]
367  - netFlux(dVf_, celli)/Vi*rDeltaT
368  - netFlux(dVfcorrectionValues, celli)/Vi*rDeltaT
369  )
370  /
371  (rDeltaT - Sp[celli]);
372 
373  alphaOvershoot =
374  pos0(alpha1New - 1)*(alpha1New - 1)
375  + neg0(alpha1New)*alpha1New;
376 
377  fluidToPassOn = alphaOvershoot*Vi;
378 
379  DebugInfo
380  << "\nNew alpha for cell " << celli << ": "
381  << alpha1New << endl;
382  }
383  }
384  }
385 
386  DebugInfo << "correctedFaces = " << correctedFaces << endl;
387 }
388 
389 
390 template<class SpType, class SuType>
391 void Foam::isoAdvection::advect(const SpType& Sp, const SuType& Su)
392 {
393  addProfilingInFunction(geometricVoF);
395 
396  if (mesh_.topoChanging())
397  {
398  setProcessorPatches();
399  }
400 
401  scalar advectionStartTime = mesh_.time().elapsedCpuTime();
402 
403  const scalar rDeltaT = 1.0/mesh_.time().deltaTValue();
404 
405  // reconstruct the interface
406  surf_->reconstruct();
407 
408  if (timeIndex_ < mesh_.time().timeIndex())
409  {
410  timeIndex_= mesh_.time().timeIndex();
411  surf_->normal().oldTime() = surf_->normal();
412  surf_->centre().oldTime() = surf_->centre();
413  }
414 
415  // Initialising dVf with upwind values
416  // i.e. phi[facei]*alpha1[upwindCell[facei]]*dt
417  dVf_ = upwind<scalar>(mesh_, phi_).flux(alpha1_)*mesh_.time().deltaT();
418 
419  // Do the isoAdvection on surface cells
420  timeIntegratedFlux();
421 
422  // Adjust alpha for mesh motion
423  if (mesh_.moving())
424  {
425  alpha1In_ *= (mesh_.Vsc0()/mesh_.Vsc());
426  }
427 
428  // Advect the free surface
429  if (porosityEnabled_)
430  {
431  // Should Su and Sp also be divided by porosity?
432  alpha1_.primitiveFieldRef() =
433  (
434  alpha1_.oldTime().primitiveField()*rDeltaT + Su.field()
435  - fvc::surfaceIntegrate(dVf_)().primitiveField()*rDeltaT
436  /porosityPtr_->primitiveField()
437  )/(rDeltaT - Sp.field());
438  }
439  else
440  {
441  alpha1_.primitiveFieldRef() =
442  (
443  alpha1_.oldTime().primitiveField()*rDeltaT
444  + Su.field()
445  - fvc::surfaceIntegrate(dVf_)().primitiveField()*rDeltaT
446  )/(rDeltaT - Sp.field());
447  }
448 
449  alpha1_.correctBoundaryConditions();
450 
451  // Adjust dVf for unbounded cells
452  limitFluxes(Sp, Su);
453 
454  scalar maxAlphaMinus1 = gMax(alpha1In_) - 1;
455  scalar minAlpha = gMin(alpha1In_);
456 
457  Info<< "isoAdvection: After conservative bounding: min(alpha) = "
458  << minAlpha << ", max(alpha) = 1 + " << maxAlphaMinus1 << endl;
459 
460  // Apply non-conservative bounding mechanisms (clipping and snapping)
461  // Note: We should be able to write out alpha before this is done!
462  applyBruteForceBounding();
463 
464  // Write surface cell set and bound cell set if required by user
465  writeSurfaceCells();
466 
467  advectionTime_ += (mesh_.time().elapsedCpuTime() - advectionStartTime);
468  DebugInfo
469  << "isoAdvection: time consumption = "
470  << label(100*advectionTime_/(mesh_.time().elapsedCpuTime() + SMALL))
471  << "%" << endl;
472 
473  alphaPhi_ = dVf_/mesh_.time().deltaT();
474 }
475 
476 
477 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
const polyBoundaryMesh & pbm
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
zeroField Su
Definition: alphaSuSp.H:1
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const labelList & patchID() const
Per boundary face label the patch index.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
dimensionedScalar neg0(const dimensionedScalar &ds)
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
#define DebugInfo
Report an information message using Foam::Info.
dimensionedScalar pos0(const dimensionedScalar &ds)
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
labelList f(nPoints)
zeroField field() const noexcept
Method name compatibility with DimensionedField.
Definition: zeroField.H:66
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void advect(const SpType &Sp, const SuType &Su)
Advect the free surface. Updates alpha field, taking into account.
#define addProfilingInFunction(Name)
Define profiling trigger with specified name and description corresponding to the compiler-defined fu...
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
List< label > labelList
A List of labels.
Definition: List.H:62
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
zeroField Sp
Definition: alphaSuSp.H:2