twoPhaseSystem.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 | Copyright (C) 2013-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "twoPhaseSystem.H"
28 #include "BlendedInterfacialModel.H"
29 #include "virtualMassModel.H"
30 #include "heatTransferModel.H"
31 #include "liftModel.H"
32 #include "wallLubricationModel.H"
33 #include "turbulentDispersionModel.H"
34 #include "fvMatrix.H"
35 #include "surfaceInterpolate.H"
36 #include "MULES.H"
37 #include "subCycle.H"
38 #include "fvcDdt.H"
39 #include "fvcDiv.H"
40 #include "fvcSnGrad.H"
41 #include "fvcFlux.H"
42 #include "fvcCurl.H"
43 #include "fvmDdt.H"
44 #include "fvmLaplacian.H"
46 
47 #include "blendingMethod.H"
48 #include "HashPtrTable.H"
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const fvMesh& mesh,
56  const dimensionedVector& g
57 )
58 :
59  IOdictionary
60  (
61  IOobject
62  (
63  "phaseProperties",
64  mesh.time().constant(),
65  mesh,
66  IOobject::MUST_READ_IF_MODIFIED,
67  IOobject::NO_WRITE
68  )
69  ),
70 
71  mesh_(mesh),
72 
73  phase1_
74  (
75  *this,
76  *this,
77  wordList(lookup("phases"))[0]
78  ),
79 
80  phase2_
81  (
82  *this,
83  *this,
84  wordList(lookup("phases"))[1]
85  ),
86 
87  phi_
88  (
89  IOobject
90  (
91  "phi",
92  mesh.time().timeName(),
93  mesh,
94  IOobject::NO_READ,
95  IOobject::AUTO_WRITE
96  ),
97  this->calcPhi()
98  ),
99 
100  dgdt_
101  (
102  IOobject
103  (
104  "dgdt",
105  mesh.time().timeName(),
106  mesh,
107  IOobject::READ_IF_PRESENT,
108  IOobject::AUTO_WRITE
109  ),
110  mesh,
111  dimensionedScalar("dgdt", dimless/dimTime, 0)
112  )
113 {
114  phase2_.volScalarField::operator=(scalar(1) - phase1_);
115 
116 
117  // Blending
118  forAllConstIter(dictionary, subDict("blending"), iter)
119  {
120  blendingMethods_.insert
121  (
122  iter().dict().dictName(),
124  (
125  iter().dict(),
126  wordList(lookup("phases"))
127  )
128  );
129  }
130 
131 
132  // Pairs
133 
134  phasePair::scalarTable sigmaTable(lookup("sigma"));
135  phasePair::dictTable aspectRatioTable(lookup("aspectRatio"));
136 
137  pair_.set
138  (
139  new phasePair
140  (
141  phase1_,
142  phase2_,
143  g,
144  sigmaTable
145  )
146  );
147 
148  pair1In2_.set
149  (
150  new orderedPhasePair
151  (
152  phase1_,
153  phase2_,
154  g,
155  sigmaTable,
156  aspectRatioTable
157  )
158  );
159 
160  pair2In1_.set
161  (
162  new orderedPhasePair
163  (
164  phase2_,
165  phase1_,
166  g,
167  sigmaTable,
168  aspectRatioTable
169  )
170  );
171 
172 
173  // Models
174 
175  drag_.set
176  (
177  new BlendedInterfacialModel<dragModel>
178  (
179  lookup("drag"),
180  (
181  blendingMethods_.found("drag")
182  ? blendingMethods_["drag"]
183  : blendingMethods_["default"]
184  ),
185  pair_,
186  pair1In2_,
187  pair2In1_,
188  false // Do not zero drag coefficent at fixed-flux BCs
189  )
190  );
191 
192  virtualMass_.set
193  (
194  new BlendedInterfacialModel<virtualMassModel>
195  (
196  lookup("virtualMass"),
197  (
198  blendingMethods_.found("virtualMass")
199  ? blendingMethods_["virtualMass"]
200  : blendingMethods_["default"]
201  ),
202  pair_,
203  pair1In2_,
204  pair2In1_
205  )
206  );
207 
208  heatTransfer_.set
209  (
210  new BlendedInterfacialModel<heatTransferModel>
211  (
212  lookup("heatTransfer"),
213  (
214  blendingMethods_.found("heatTransfer")
215  ? blendingMethods_["heatTransfer"]
216  : blendingMethods_["default"]
217  ),
218  pair_,
219  pair1In2_,
220  pair2In1_
221  )
222  );
223 
224  lift_.set
225  (
226  new BlendedInterfacialModel<liftModel>
227  (
228  lookup("lift"),
229  (
230  blendingMethods_.found("lift")
231  ? blendingMethods_["lift"]
232  : blendingMethods_["default"]
233  ),
234  pair_,
235  pair1In2_,
236  pair2In1_
237  )
238  );
239 
240  wallLubrication_.set
241  (
242  new BlendedInterfacialModel<wallLubricationModel>
243  (
244  lookup("wallLubrication"),
245  (
246  blendingMethods_.found("wallLubrication")
247  ? blendingMethods_["wallLubrication"]
248  : blendingMethods_["default"]
249  ),
250  pair_,
251  pair1In2_,
252  pair2In1_
253  )
254  );
255 
256  turbulentDispersion_.set
257  (
258  new BlendedInterfacialModel<turbulentDispersionModel>
259  (
260  lookup("turbulentDispersion"),
261  (
262  blendingMethods_.found("turbulentDispersion")
263  ? blendingMethods_["turbulentDispersion"]
264  : blendingMethods_["default"]
265  ),
266  pair_,
267  pair1In2_,
268  pair2In1_
269  )
270  );
271 }
272 
273 
274 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
275 
277 {}
278 
279 
280 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
281 
283 {
284  return phase1_*phase1_.thermo().rho() + phase2_*phase2_.thermo().rho();
285 }
286 
287 
289 {
290  return phase1_*phase1_.U() + phase2_*phase2_.U();
291 }
292 
293 
295 {
296  return
297  fvc::interpolate(phase1_)*phase1_.phi()
298  + fvc::interpolate(phase2_)*phase2_.phi();
299 }
300 
301 
303 {
304  return drag_->K();
305 }
306 
307 
309 {
310  return drag_->Kf();
311 }
312 
313 
315 {
316  return virtualMass_->K();
317 }
318 
319 
321 {
322  return virtualMass_->Kf();
323 }
324 
325 
327 {
328  return heatTransfer_->K();
329 }
330 
331 
333 {
334  return lift_->F<vector>() + wallLubrication_->F<vector>();
335 }
336 
337 
339 {
340  return lift_->Ff() + wallLubrication_->Ff();
341 }
342 
343 
345 {
346  return turbulentDispersion_->D();
347 }
348 
349 
351 {
352  const Time& runTime = mesh_.time();
353 
354  volScalarField& alpha1 = phase1_;
355  volScalarField& alpha2 = phase2_;
356 
357  const surfaceScalarField& phi1 = phase1_.phi();
358  const surfaceScalarField& phi2 = phase2_.phi();
359 
360  const dictionary& alphaControls = mesh_.solverDict
361  (
362  alpha1.name()
363  );
364 
365  label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")));
366  label nAlphaCorr(readLabel(alphaControls.lookup("nAlphaCorr")));
367 
368  word alphaScheme("div(phi," + alpha1.name() + ')');
369  word alpharScheme("div(phir," + alpha1.name() + ')');
370 
371  alpha1.correctBoundaryConditions();
372 
373 
374  surfaceScalarField phic("phic", phi_);
375  surfaceScalarField phir("phir", phi1 - phi2);
376 
377  tmp<surfaceScalarField> alpha1alpha2f;
378 
379  if (pPrimeByA_.valid())
380  {
381  alpha1alpha2f =
382  fvc::interpolate(max(alpha1, scalar(0)))
383  *fvc::interpolate(max(alpha2, scalar(0)));
384 
385  surfaceScalarField phiP
386  (
387  pPrimeByA_()*fvc::snGrad(alpha1, "bounded")*mesh_.magSf()
388  );
389 
390  phir += phiP;
391  }
392 
393  for (int acorr=0; acorr<nAlphaCorr; acorr++)
394  {
396  (
397  IOobject
398  (
399  "Sp",
400  runTime.timeName(),
401  mesh_
402  ),
403  mesh_,
404  dimensionedScalar("Sp", dgdt_.dimensions(), 0.0)
405  );
406 
408  (
409  IOobject
410  (
411  "Su",
412  runTime.timeName(),
413  mesh_
414  ),
415  // Divergence term is handled explicitly to be
416  // consistent with the explicit transport solution
417  fvc::div(phi_)*min(alpha1, scalar(1))
418  );
419 
420  forAll(dgdt_, celli)
421  {
422  if (dgdt_[celli] > 0.0)
423  {
424  Sp[celli] -= dgdt_[celli]/max(1.0 - alpha1[celli], 1e-4);
425  Su[celli] += dgdt_[celli]/max(1.0 - alpha1[celli], 1e-4);
426  }
427  else if (dgdt_[celli] < 0.0)
428  {
429  Sp[celli] += dgdt_[celli]/max(alpha1[celli], 1e-4);
430  }
431  }
432 
433  surfaceScalarField alphaPhic1
434  (
435  fvc::flux
436  (
437  phic,
438  alpha1,
439  alphaScheme
440  )
441  + fvc::flux
442  (
443  -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
444  alpha1,
446  )
447  );
448 
449  // Ensure that the flux at inflow BCs is preserved
450  forAll(alphaPhic1.boundaryField(), patchi)
451  {
452  fvsPatchScalarField& alphaPhic1p =
453  alphaPhic1.boundaryField()[patchi];
454 
455  if (!alphaPhic1p.coupled())
456  {
457  const scalarField& phi1p = phi1.boundaryField()[patchi];
458  const scalarField& alpha1p = alpha1.boundaryField()[patchi];
459 
460  forAll(alphaPhic1p, facei)
461  {
462  if (phi1p[facei] < 0)
463  {
464  alphaPhic1p[facei] = alpha1p[facei]*phi1p[facei];
465  }
466  }
467  }
468  }
469 
470  if (nAlphaSubCycles > 1)
471  {
472  for
473  (
474  subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
475  !(++alphaSubCycle).end();
476  )
477  {
478  surfaceScalarField alphaPhic10(alphaPhic1);
479 
481  (
482  geometricOneField(),
483  alpha1,
484  phi_,
485  alphaPhic10,
486  (alphaSubCycle.index()*Sp)(),
487  (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
488  phase1_.alphaMax(),
489  0
490  );
491 
492  if (alphaSubCycle.index() == 1)
493  {
494  phase1_.alphaPhi() = alphaPhic10;
495  }
496  else
497  {
498  phase1_.alphaPhi() += alphaPhic10;
499  }
500  }
501 
502  phase1_.alphaPhi() /= nAlphaSubCycles;
503  }
504  else
505  {
507  (
508  geometricOneField(),
509  alpha1,
510  phi_,
511  alphaPhic1,
512  Sp,
513  Su,
514  phase1_.alphaMax(),
515  0
516  );
517 
518  phase1_.alphaPhi() = alphaPhic1;
519  }
520 
521  if (pPrimeByA_.valid())
522  {
523  fvScalarMatrix alpha1Eqn
524  (
526  - fvm::laplacian(alpha1alpha2f*pPrimeByA_(), alpha1, "bounded")
527  );
528 
529  alpha1Eqn.relax();
530  alpha1Eqn.solve();
531 
532  phase1_.alphaPhi() += alpha1Eqn.flux();
533  }
534 
535  phase1_.alphaRhoPhi() =
536  fvc::interpolate(phase1_.rho())*phase1_.alphaPhi();
537 
538  phase2_.alphaPhi() = phi_ - phase1_.alphaPhi();
539  alpha2 = scalar(1) - alpha1;
540  phase2_.alphaRhoPhi() =
541  fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();
542 
543  Info<< alpha1.name() << " volume fraction = "
544  << alpha1.weightedAverage(mesh_.V()).value()
545  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
546  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
547  << endl;
548  }
549 }
550 
551 
553 {
554  phase1_.correct();
555  phase2_.correct();
556 }
557 
558 
560 {
561  phase1_.turbulence().correct();
562  phase2_.turbulence().correct();
563 }
564 
565 
567 {
568  if (regIOobject::read())
569  {
570  bool readOK = true;
571 
572  readOK &= phase1_.read(*this);
573  readOK &= phase2_.read(*this);
574 
575  // models ...
576 
577  return readOK;
578  }
579  else
580  {
581  return false;
582  }
583 }
584 
585 
586 const Foam::dragModel& Foam::twoPhaseSystem::drag(const phaseModel& phase) const
587 {
588  return drag_->phaseModel(phase);
589 }
590 
591 
593 Foam::twoPhaseSystem::virtualMass(const phaseModel& phase) const
594 {
595  return virtualMass_->phaseModel(phase);
596 }
597 
598 
600 {
601  return pair_->sigma();
602 }
603 
604 
605 // ************************************************************************* //
Foam::virtualMassModel
Definition: virtualMassModel.H:52
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::MULES::explicitSolve
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Definition: MULESTemplates.C:38
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
phi1
surfaceScalarField & phi1
Definition: createFields.H:19
Foam::twoPhaseSystem::rho
tmp< volScalarField > rho() const
Return the mixture density.
Foam::blendingMethod::New
static autoPtr< blendingMethod > New(const dictionary &dict, const wordList &phaseNames)
Foam::twoPhaseSystem::Vm
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
subCycle.H
Foam::fvsPatchScalarField
fvsPatchField< scalar > fvsPatchScalarField
Definition: fvsPatchFieldsFwd.H:43
fvcSnGrad.H
Calculate the snGrad of the given volField.
Foam::twoPhaseSystem::correct
void correct()
Correct two-phase properties other than turbulence.
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::regIOobject::read
virtual bool read()
Read object.
Definition: regIOobjectRead.C:171
fvcDiv.H
Calculate the divergence of the given field.
Foam::twoPhaseSystem::Kh
tmp< volScalarField > Kh() const
Return the heat transfer coefficient.
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::twoPhaseSystem::sigma
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
phir
surfaceScalarField phir(phic *interface.nHatf())
phi2
surfaceScalarField & phi2
Definition: createFields.H:24
Foam::twoPhaseSystem::Kd
tmp< volScalarField > Kd() const
Return the drag coefficient.
PhaseCompressibleTurbulenceModel.H
constant
Constant dispersed-phase particle diameter model.
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
fvcCurl.H
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
Foam::twoPhaseSystem::calcPhi
tmp< surfaceScalarField > calcPhi() const
Return the mixture flux.
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
dictName
const word dictName("particleTrackDict")
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
nAlphaSubCycles
label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")))
Foam::Info
messageStream Info
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
alpha2
alpha2
Definition: alphaEqn.H:112
dict
dictionary dict
Definition: searchingEngine.H:14
fixedValueFvsPatchFields.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::twoPhaseSystem::twoPhaseSystem
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
surfaceInterpolate.H
Surface Interpolation.
Foam::phasePair::scalarTable
HashTable< scalar, phasePairKey, phasePairKey::hash > scalarTable
Scalar hash table.
Definition: phasePair.H:63
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
phic
phic
Definition: alphaEqnsSubCycle.H:5
alpha1
volScalarField & alpha1
Definition: createFields.H:15
Foam::twoPhaseSystem::Vmf
tmp< surfaceScalarField > Vmf() const
Return the face virtual mass coefficient.
Foam::twoPhaseSystem::virtualMass
const virtualMassModel & virtualMass(const phaseModel &phase) const
Return the virtual mass model for the given phase.
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
alpharScheme
word alpharScheme("div(phirb,alpha)")
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::twoPhaseSystem::Kdf
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
Foam::twoPhaseSystem::D
tmp< volScalarField > D() const
Return the turbulent diffusivity.
Foam::phasePair::dictTable
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition: phasePair.H:59
Foam::twoPhaseSystem::read
bool read()
Read base phaseProperties dictionary.
fvcFlux.H
Calculate the face-flux of the given field.
nAlphaCorr
label nAlphaCorr(readLabel(alphaControls.lookup("nAlphaCorr")))
patchi
label patchi
Definition: getPatchFieldScalar.H:1
HashPtrTable.H
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::twoPhaseSystem::correctTurbulence
void correctTurbulence()
Correct two-phase turbulence.
fvcDdt.H
Calculate the first temporal derivative.
Foam::GeometricField::DimensionedInternalField
DimensionedField< Type, GeoMesh > DimensionedInternalField
Definition: GeometricField.H:100
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
Foam::dragModel
Definition: dragModel.H:50
Foam::twoPhaseSystem::F
tmp< volVectorField > F() const
Return the combined force (lift + wall-lubrication)
Foam::twoPhaseSystem::solve
virtual void solve()
Solve for the phase fractions.
Foam::twoPhaseSystem::U
tmp< volVectorField > U() const
Return the mixture velocity.
Foam::fvc::flux
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcFlux.C:45
Foam::fvc::Sp
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
MULES.H
MULES: Multidimensional universal limiter for explicit solution.
Foam::twoPhaseSystem::drag
const dragModel & drag(const phaseModel &phase) const
Return the drag model for the given phase.
Foam::twoPhaseSystem::~twoPhaseSystem
virtual ~twoPhaseSystem()
Destructor.
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fvmDdt.H
Calulate the matrix for the first temporal derivative.
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
alphaControls
const dictionary & alphaControls
Definition: alphaControls.H:1
Foam::twoPhaseSystem::Ff
tmp< surfaceScalarField > Ff() const
Return the combined face-force (lift + wall-lubrication)