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 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "isoAdvection.H"
31 #include "fvcSurfaceIntegrate.H"
32 #include "upwind.H"
33 
34 // ************************************************************************* //
35 
36 template<typename Type>
37 Type Foam::isoAdvection::faceValue
38 (
39  const GeometricField<Type, fvsPatchField, surfaceMesh>& f,
40  const label facei
41 ) const
42 {
43  if (mesh_.isInternalFace(facei))
44  {
45  return f.primitiveField()[facei];
46  }
47  else
48  {
49  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
50 
51  // Boundary face. Find out which face of which patch
52  const label patchi = pbm.patchID()[facei - mesh_.nInternalFaces()];
53 
54  if (patchi < 0 || patchi >= pbm.size())
55  {
57  << "Cannot find patch for face " << facei
58  << abort(FatalError);
59  }
60 
61  // Handle empty patches
62  const polyPatch& pp = pbm[patchi];
63  if (isA<emptyPolyPatch>(pp) || pp.empty())
64  {
65  return pTraits<Type>::zero;
66  }
67 
68  const label patchFacei = pp.whichFace(facei);
69  return f.boundaryField()[patchi][patchFacei];
70  }
71 }
72 
73 
74 template<typename Type>
75 void Foam::isoAdvection::setFaceValue
76 (
77  GeometricField<Type, fvsPatchField, surfaceMesh>& f,
78  const label facei,
79  const Type& value
80 ) const
81 {
82  if (mesh_.isInternalFace(facei))
83  {
84  f.primitiveFieldRef()[facei] = value;
85  }
86  else
87  {
88  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
89 
90  // Boundary face. Find out which face of which patch
91  const label patchi = pbm.patchID()[facei - mesh_.nInternalFaces()];
92 
93  if (patchi < 0 || patchi >= pbm.size())
94  {
96  << "Cannot find patch for face " << facei
97  << abort(FatalError);
98  }
99 
100  // Handle empty patches
101  const polyPatch& pp = pbm[patchi];
102  if (isA<emptyPolyPatch>(pp) || pp.empty())
103  {
104  return;
105  }
106 
107  const label patchFacei = pp.whichFace(facei);
108 
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 {
122 
123  const scalar aTol = 1.0e-12; // Note: tolerances
124  scalar maxAlphaMinus1 = gMax(alpha1_) - 1; // max(alphaNew - 1);
125  scalar minAlpha = gMin(alpha1_); // min(alphaNew);
126  const label nOvershoots = 20; // sum(pos0(alphaNew - 1 - aTol));
127 
128  const labelList& owner = mesh_.faceOwner();
129  const labelList& neighbour = mesh_.faceNeighbour();
130 
131  DebugInfo
132  << "isoAdvection: Before conservative bounding: min(alpha) = "
133  << minAlpha << ", max(alpha) = 1 + " << maxAlphaMinus1 << endl;
134 
135  surfaceScalarField dVfcorrectionValues("dVfcorrectionValues", dVf_*0.0);
136 
137 
138  // Loop number of bounding steps
139  for (label n = 0; n < nAlphaBounds_; n++)
140  {
141  if (maxAlphaMinus1 > aTol || minAlpha < -aTol) // Note: tolerances
142  {
143  DebugInfo << "boundAlpha... " << endl;
144 
145  DynamicList<label> correctedFaces(3*nOvershoots);
146  dVfcorrectionValues = dimensionedScalar(dimVolume, Zero);
147  boundFlux(alpha1In_, dVfcorrectionValues, correctedFaces,Sp,Su);
148 
149  correctedFaces.append
150  (
151  syncProcPatches(dVfcorrectionValues, phi_,true)
152  );
153 
154  labelHashSet alreadyUpdated;
155  forAll(correctedFaces, fi)
156  {
157  label facei = correctedFaces[fi];
158  if (alreadyUpdated.insert(facei))
159  {
160  checkIfOnProcPatch(facei);
161  const label own = owner[facei];
162 
163  alpha1_[own] +=
164  -faceValue(dVfcorrectionValues, facei)/mesh_.V()[own];
165  if (mesh_.isInternalFace(facei))
166  {
167  const label nei = neighbour[facei];
168  alpha1_[nei] -=
169  -faceValue(dVfcorrectionValues, facei)/mesh_.V()[nei];
170  }
171 
172  // Change to treat boundaries consistently
173  scalar corrVf =
174  faceValue(dVf_, facei)
175  + faceValue(dVfcorrectionValues, facei);
176 
177  setFaceValue(dVf_, facei, corrVf);
178  }
179  }
180  syncProcPatches(dVf_, phi_);
181  }
182  else
183  {
184  break;
185  }
186 
187  maxAlphaMinus1 = gMax(alpha1_) - 1; // max(alphaNew - 1);
188  minAlpha = gMin(alpha1_); // min(alphaNew);
189 
190  if (debug)
191  {
192  // Check if still unbounded
193  //scalarField alphaNew(alpha1In_ - fvc::surfaceIntegrate(dVf_)());
194  label maxAlphaMinus1 = max(alpha1_.primitiveField() - 1);
195  scalar minAlpha = min(alpha1_.primitiveField());
196  label nUndershoots = sum(neg0(alpha1_.primitiveField() + aTol));
197  label nOvershoots = sum(pos0(alpha1_.primitiveField() - 1 - aTol));
198 
199  Info<< "After bounding number " << n + 1 << " of time "
200  << mesh_.time().value() << ":" << nl
201  << "nOvershoots = " << nOvershoots << " with max(alpha1_-1) = "
202  << maxAlphaMinus1 << " and nUndershoots = " << nUndershoots
203  << " with min(alpha1_) = " << minAlpha << endl;
204  }
205  }
206 
207  alpha1_.correctBoundaryConditions();
208 
209 }
210 
211 
212 template<class SpType, class SuType>
213 void Foam::isoAdvection::boundFlux
214 (
215  const scalarField& alpha1,
216  surfaceScalarField& dVfcorrectionValues,
217  DynamicList<label>& correctedFaces,
218  const SpType& Sp,
219  const SuType& Su
220 )
221 {
223  scalar rDeltaT = 1/mesh_.time().deltaTValue();
224 
225  correctedFaces.clear();
226  scalar aTol = 10*SMALL; // Note: tolerances
227 
228  const scalarField& meshV = mesh_.cellVolumes();
229  const scalar dt = mesh_.time().deltaTValue();
230 
231  DynamicList<label> downwindFaces(10);
232  DynamicList<label> facesToPassFluidThrough(downwindFaces.size());
233  DynamicList<scalar> dVfmax(downwindFaces.size());
234  DynamicList<scalar> phi(downwindFaces.size());
235  const volScalarField& alphaOld = alpha1_.oldTime();
236 
237  // Loop through alpha cell centred field
238  forAll(alpha1, celli)
239  {
240  if (alpha1[celli] < -aTol || alpha1[celli] > 1 + aTol)
241  {
242  const scalar Vi = meshV[celli];
243  scalar alphaOvershoot = pos0(alpha1[celli]-1)*(alpha1[celli]-1)
244  + neg0(alpha1[celli])*alpha1[celli];
245  scalar fluidToPassOn = alphaOvershoot*Vi;
246  label nFacesToPassFluidThrough = 1;
247 
248  bool firstLoop = true;
249 
250  // First try to pass surplus fluid on to neighbour cells that are
251  // not filled and to which dVf < phi*dt
252  while (mag(alphaOvershoot) > aTol && nFacesToPassFluidThrough > 0)
253  {
254  DebugInfo
255  << "\n\nBounding cell " << celli
256  << " with alpha overshooting " << alphaOvershoot
257  << endl;
258 
259  facesToPassFluidThrough.clear();
260  dVfmax.clear();
261  phi.clear();
262 
263  // Find potential neighbour cells to pass surplus phase to
264  setDownwindFaces(celli, downwindFaces);
265 
266  scalar dVftot = 0;
267  nFacesToPassFluidThrough = 0;
268 
269  forAll(downwindFaces, fi)
270  {
271  const label facei = downwindFaces[fi];
272  const scalar phif = faceValue(phi_, facei);
273 
274  const scalar dVff =
275  faceValue(dVf_, facei)
276  + faceValue(dVfcorrectionValues, facei);
277 
278  const scalar maxExtraFaceFluidTrans =
279  mag(pos0(fluidToPassOn)*phif*dt - dVff);
280 
281  // dVf has same sign as phi and so if phi>0 we have
282  // mag(phi_[facei]*dt) - mag(dVf[facei]) = phi_[facei]*dt
283  // - dVf[facei]
284  // If phi < 0 we have mag(phi_[facei]*dt) -
285  // mag(dVf[facei]) = -phi_[facei]*dt - (-dVf[facei]) > 0
286  // since mag(dVf) < phi*dt
287  DebugInfo
288  << "downwindFace " << facei
289  << " has maxExtraFaceFluidTrans = "
290  << maxExtraFaceFluidTrans << endl;
291 
292  if (maxExtraFaceFluidTrans/Vi > aTol)
293  {
294 // if (maxExtraFaceFluidTrans/Vi > aTol &&
295 // mag(dVfIn[facei])/Vi > aTol) //Last condition may be
296 // important because without this we will flux through uncut
297 // downwind faces
298  facesToPassFluidThrough.append(facei);
299  phi.append(phif);
300  dVfmax.append(maxExtraFaceFluidTrans);
301  dVftot += mag(phif*dt);
302  }
303  }
304 
305  DebugInfo
306  << "\nfacesToPassFluidThrough: "
307  << facesToPassFluidThrough << ", dVftot = "
308  << dVftot << " m3 corresponding to dalpha = "
309  << dVftot/Vi << endl;
310 
311  forAll(facesToPassFluidThrough, fi)
312  {
313  const label facei = facesToPassFluidThrough[fi];
314  scalar fluidToPassThroughFace =
315  mag(fluidToPassOn)*mag(phi[fi]*dt)/dVftot;
316 
317  nFacesToPassFluidThrough +=
318  pos0(dVfmax[fi] - fluidToPassThroughFace);
319 
320  fluidToPassThroughFace =
321  min(fluidToPassThroughFace, dVfmax[fi]);
322 
323  scalar dVff = faceValue(dVfcorrectionValues, facei);
324 
325  dVff +=
326  sign(phi[fi])*fluidToPassThroughFace*sign(fluidToPassOn);
327 
328  setFaceValue(dVfcorrectionValues, facei, dVff);
329 
330  if (firstLoop)
331  {
332  checkIfOnProcPatch(facei);
333  correctedFaces.append(facei);
334  }
335  }
336 
337  firstLoop = false;
338 
339  scalar alpha1New =
340  (
341  alphaOld[celli]*rDeltaT + Su[celli]
342  - netFlux(dVf_, celli)/Vi*rDeltaT
343  - netFlux(dVfcorrectionValues, celli)/Vi*rDeltaT
344  )
345  /
346  (rDeltaT - Sp[celli]);
347 
348  alphaOvershoot =
349  pos0(alpha1New-1)*(alpha1New-1)
350  + neg0(alpha1New)*alpha1New;
351 
352  fluidToPassOn = alphaOvershoot*Vi;
353 
354  DebugInfo
355  << "\nNew alpha for cell " << celli << ": "
356  << alpha1New << endl;
357  }
358  }
359  }
360 
361  DebugInfo << "correctedFaces = " << correctedFaces << endl;
362 }
363 
364 
365 template<class SpType, class SuType>
366 void Foam::isoAdvection::advect(const SpType& Sp, const SuType& Su)
367 {
369 
370  scalar advectionStartTime = mesh_.time().elapsedCpuTime();
371 
372  scalar rDeltaT = 1/mesh_.time().deltaTValue();
373 
374  // reconstruct the interface
375  surf_->reconstruct();
376 
377  // Initialising dVf with upwind values
378  // i.e. phi[facei]*alpha1[upwindCell[facei]]*dt
379  dVf_ = upwind<scalar>(mesh_, phi_).flux(alpha1_)*mesh_.time().deltaT();
380 
381  // Do the isoAdvection on surface cells
382  timeIntegratedFlux();
383 
384  // Adjust alpha for mesh motion
385  if (mesh_.moving())
386  {
387  alpha1In_ *= (mesh_.Vsc0()/mesh_.Vsc());
388  }
389 
390  // Advect the free surface
391  alpha1_.primitiveFieldRef() =
392  (
393  alpha1_.oldTime().primitiveField()*rDeltaT
394  + Su.field()
395  - fvc::surfaceIntegrate(dVf_)().primitiveField()*rDeltaT
396  )/(rDeltaT - Sp.field());
397 
398  alpha1_.correctBoundaryConditions();
399 
400  // Adjust dVf for unbounded cells
401  limitFluxes
402  (
403  Sp,
404  Su
405  );
406 
407  scalar maxAlphaMinus1 = gMax(alpha1In_) - 1;
408  scalar minAlpha = gMin(alpha1In_);
409 
410  Info<< "isoAdvection: After conservative bounding: min(alpha) = "
411  << minAlpha << ", max(alpha) = 1 + " << maxAlphaMinus1 << endl;
412 
413  // Apply non-conservative bounding mechanisms (clipping and snapping)
414  // Note: We should be able to write out alpha before this is done!
415  applyBruteForceBounding();
416 
417  // Write surface cell set and bound cell set if required by user
418  writeSurfaceCells();
419 
420  advectionTime_ += (mesh_.time().elapsedCpuTime() - advectionStartTime);
421  DebugInfo
422  << "isoAdvection: time consumption = "
423  << label(100*advectionTime_/(mesh_.time().elapsedCpuTime() + SMALL))
424  << "%" << endl;
425 
426  alphaPhi_ = dVf_/mesh_.time().deltaT();
427 }
428 
429 
430 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::upwind
Upwind differencing scheme class.
Definition: upwind.H:56
Foam::DynamicList< label >
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::neg0
dimensionedScalar neg0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:210
Sp
zeroField Sp
Definition: alphaSuSp.H:2
Foam::GeometricField::oldTime
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Definition: GeometricField.C:850
Foam::isoAdvection::advect
void advect(const SpType &Sp, const SuType &Su)
Advect the free surface. Updates alpha field, taking into account.
Definition: isoAdvectionTemplates.C:366
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Su
zeroField Su
Definition: alphaSuSp.H:1
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::polyBoundaryMesh::patchID
const labelList & patchID() const
Per boundary face label the patch index.
Definition: polyBoundaryMesh.C:452
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:474
Foam::DynamicList::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:348
upwind.H
Foam::zeroField::field
zeroField field() const noexcept
Definition: zeroField.H:67
Foam::limitedSurfaceInterpolationScheme::flux
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Definition: limitedSurfaceInterpolationScheme.C:195
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
fvcSurfaceIntegrate.H
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
isoAdvection.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:359
Foam::nl
constexpr char nl
Definition: Ostream.H:385
f
labelList f(nPoints)
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::fvc::surfaceIntegrate
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcSurfaceIntegrate.C:46
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:102
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:409
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592