CMULESTemplates.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-2017 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 "CMULES.H"
27 #include "fvcSurfaceIntegrate.H"
28 #include "localEulerDdtScheme.H"
29 #include "slicedSurfaceFields.H"
30 #include "wedgeFvPatch.H"
31 #include "syncTools.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 template<class RdeltaTType, class RhoType, class SpType, class SuType>
37 (
38  const RdeltaTType& rDeltaT,
39  const RhoType& rho,
40  volScalarField& psi,
41  const surfaceScalarField& phi,
42  const surfaceScalarField& phiCorr,
43  const SpType& Sp,
44  const SuType& Su
45 )
46 {
47  Info<< "MULES: Correcting " << psi.name() << endl;
48 
49  const fvMesh& mesh = psi.mesh();
50 
51  scalarField psiIf(psi.size(), 0);
52  fvc::surfaceIntegrate(psiIf, phiCorr);
53 
54  if (mesh.moving())
55  {
56  psi.primitiveFieldRef() =
57  (
58  rho.field()*psi.primitiveField()*rDeltaT
59  + Su.field()
60  - psiIf
61  )/(rho.field()*rDeltaT - Sp.field());
62  }
63  else
64  {
65  psi.primitiveFieldRef() =
66  (
67  rho.field()*psi.primitiveField()*rDeltaT
68  + Su.field()
69  - psiIf
70  )/(rho.field()*rDeltaT - Sp.field());
71  }
72 
74 }
75 
76 
77 template<class RhoType, class SpType, class SuType>
79 (
80  const RhoType& rho,
81  volScalarField& psi,
82  const surfaceScalarField& phi,
83  surfaceScalarField& phiCorr,
84  const SpType& Sp,
85  const SuType& Su,
86  const scalar psiMax,
87  const scalar psiMin
88 )
89 {
90  const fvMesh& mesh = psi.mesh();
91 
92  if (fv::localEulerDdt::enabled(mesh))
93  {
94  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
95 
96  limitCorr
97  (
98  rDeltaT,
99  rho,
100  psi,
101  phi,
102  phiCorr,
103  Sp,
104  Su,
105  psiMax,
106  psiMin
107  );
108  correct(rDeltaT, rho, psi, phi, phiCorr, Sp, Su);
109  }
110  else
111  {
112  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
113 
114  limitCorr
115  (
116  rDeltaT,
117  rho,
118  psi,
119  phi,
120  phiCorr,
121  Sp,
122  Su,
123  psiMax,
124  psiMin
125  );
126 
127  correct(rDeltaT, rho, psi, phi, phiCorr, Sp, Su);
128  }
129 }
130 
131 
132 template<class RdeltaTType, class RhoType, class SpType, class SuType>
134 (
135  scalarField& allLambda,
136  const RdeltaTType& rDeltaT,
137  const RhoType& rho,
138  const volScalarField& psi,
139  const surfaceScalarField& phi,
140  const surfaceScalarField& phiCorr,
141  const SpType& Sp,
142  const SuType& Su,
143  const scalar psiMax,
144  const scalar psiMin
145 )
146 {
147  const scalarField& psiIf = psi;
148  const volScalarField::Boundary& psiBf = psi.boundaryField();
149 
150  const fvMesh& mesh = psi.mesh();
151 
152  const dictionary& MULEScontrols = mesh.solverDict(psi.name());
153 
154  const label nLimiterIter
155  (
156  readLabel(MULEScontrols.lookup("nLimiterIter"))
157  );
158 
159  const scalar smoothLimiter
160  (
161  MULEScontrols.lookupOrDefault<scalar>("smoothLimiter", 0)
162  );
163 
164  const scalar extremaCoeff
165  (
166  MULEScontrols.lookupOrDefault<scalar>("extremaCoeff", 0)
167  );
168 
169  const labelUList& owner = mesh.owner();
170  const labelUList& neighb = mesh.neighbour();
171  tmp<volScalarField::Internal> tVsc = mesh.Vsc();
172  const scalarField& V = tVsc();
173 
174  const surfaceScalarField::Boundary& phiBf =
175  phi.boundaryField();
176 
177  const scalarField& phiCorrIf = phiCorr;
178  const surfaceScalarField::Boundary& phiCorrBf =
179  phiCorr.boundaryField();
180 
182  (
183  IOobject
184  (
185  "lambda",
186  mesh.time().timeName(),
187  mesh,
188  IOobject::NO_READ,
189  IOobject::NO_WRITE,
190  false
191  ),
192  mesh,
193  dimless,
194  allLambda,
195  false // Use slices for the couples
196  );
197 
198  scalarField& lambdaIf = lambda;
199  surfaceScalarField::Boundary& lambdaBf =
200  lambda.boundaryFieldRef();
201 
202  scalarField psiMaxn(psiIf.size(), psiMin);
203  scalarField psiMinn(psiIf.size(), psiMax);
204 
205  scalarField sumPhip(psiIf.size(), 0.0);
206  scalarField mSumPhim(psiIf.size(), 0.0);
207 
208  forAll(phiCorrIf, facei)
209  {
210  const label own = owner[facei];
211  const label nei = neighb[facei];
212 
213  psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
214  psiMinn[own] = min(psiMinn[own], psiIf[nei]);
215 
216  psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
217  psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
218 
219  const scalar phiCorrf = phiCorrIf[facei];
220 
221  if (phiCorrf > 0)
222  {
223  sumPhip[own] += phiCorrf;
224  mSumPhim[nei] += phiCorrf;
225  }
226  else
227  {
228  mSumPhim[own] -= phiCorrf;
229  sumPhip[nei] -= phiCorrf;
230  }
231  }
232 
233  forAll(phiCorrBf, patchi)
234  {
235  const fvPatchScalarField& psiPf = psiBf[patchi];
236  const scalarField& phiCorrPf = phiCorrBf[patchi];
237 
238  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
239 
240  if (psiPf.coupled())
241  {
242  const scalarField psiPNf(psiPf.patchNeighbourField());
243 
244  forAll(phiCorrPf, pFacei)
245  {
246  label pfCelli = pFaceCells[pFacei];
247 
248  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
249  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
250  }
251  }
252  else if (psiPf.fixesValue())
253  {
254  forAll(phiCorrPf, pFacei)
255  {
256  const label pfCelli = pFaceCells[pFacei];
257 
258  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
259  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
260  }
261  }
262  else
263  {
264  forAll(phiCorrPf, pFacei)
265  {
266  const label pfCelli = pFaceCells[pFacei];
267 
268  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiMax);
269  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiMin);
270  }
271  }
272 
273  forAll(phiCorrPf, pFacei)
274  {
275  const label pfCelli = pFaceCells[pFacei];
276 
277  const scalar phiCorrf = phiCorrPf[pFacei];
278 
279  if (phiCorrf > 0)
280  {
281  sumPhip[pfCelli] += phiCorrf;
282  }
283  else
284  {
285  mSumPhim[pfCelli] -= phiCorrf;
286  }
287  }
288  }
289 
290  psiMaxn = min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
291  psiMinn = max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
292 
293  if (smoothLimiter > SMALL)
294  {
295  psiMaxn =
296  min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
297  psiMinn =
298  max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
299  }
300 
301  psiMaxn =
302  V
303  *(
304  (rho.field()*rDeltaT - Sp.field())*psiMaxn
305  - Su.field()
306  - rho.field()*psi.primitiveField()*rDeltaT
307  );
308 
309  psiMinn =
310  V
311  *(
312  Su.field()
313  - (rho.field()*rDeltaT - Sp.field())*psiMinn
314  + rho.field()*psi.primitiveField()*rDeltaT
315  );
316 
317  scalarField sumlPhip(psiIf.size());
318  scalarField mSumlPhim(psiIf.size());
319 
320  for (int j=0; j<nLimiterIter; j++)
321  {
322  sumlPhip = 0;
323  mSumlPhim = 0;
324 
325  forAll(lambdaIf, facei)
326  {
327  const label own = owner[facei];
328  const label nei = neighb[facei];
329 
330  const scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
331 
332  if (lambdaPhiCorrf > 0)
333  {
334  sumlPhip[own] += lambdaPhiCorrf;
335  mSumlPhim[nei] += lambdaPhiCorrf;
336  }
337  else
338  {
339  mSumlPhim[own] -= lambdaPhiCorrf;
340  sumlPhip[nei] -= lambdaPhiCorrf;
341  }
342  }
343 
344  forAll(lambdaBf, patchi)
345  {
346  scalarField& lambdaPf = lambdaBf[patchi];
347  const scalarField& phiCorrfPf = phiCorrBf[patchi];
348 
349  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
350 
351  forAll(lambdaPf, pFacei)
352  {
353  label pfCelli = pFaceCells[pFacei];
354 
355  scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
356 
357  if (lambdaPhiCorrf > 0)
358  {
359  sumlPhip[pfCelli] += lambdaPhiCorrf;
360  }
361  else
362  {
363  mSumlPhim[pfCelli] -= lambdaPhiCorrf;
364  }
365  }
366  }
367 
368  forAll(sumlPhip, celli)
369  {
370  sumlPhip[celli] =
371  max(min
372  (
373  (sumlPhip[celli] + psiMaxn[celli])
374  /(mSumPhim[celli] + ROOTVSMALL),
375  1.0), 0.0
376  );
377 
378  mSumlPhim[celli] =
379  max(min
380  (
381  (mSumlPhim[celli] + psiMinn[celli])
382  /(sumPhip[celli] + ROOTVSMALL),
383  1.0), 0.0
384  );
385  }
386 
387  const scalarField& lambdam = sumlPhip;
388  const scalarField& lambdap = mSumlPhim;
389 
390  forAll(lambdaIf, facei)
391  {
392  if (phiCorrIf[facei] > 0)
393  {
394  lambdaIf[facei] = min
395  (
396  lambdaIf[facei],
397  min(lambdap[owner[facei]], lambdam[neighb[facei]])
398  );
399  }
400  else
401  {
402  lambdaIf[facei] = min
403  (
404  lambdaIf[facei],
405  min(lambdam[owner[facei]], lambdap[neighb[facei]])
406  );
407  }
408  }
409 
410 
411  forAll(lambdaBf, patchi)
412  {
413  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
414  const scalarField& phiCorrfPf = phiCorrBf[patchi];
415  const fvPatchScalarField& psiPf = psiBf[patchi];
416 
417  if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
418  {
419  lambdaPf = 0;
420  }
421  else if (psiPf.coupled())
422  {
423  const labelList& pFaceCells =
424  mesh.boundary()[patchi].faceCells();
425 
426  forAll(lambdaPf, pFacei)
427  {
428  const label pfCelli = pFaceCells[pFacei];
429 
430  if (phiCorrfPf[pFacei] > 0)
431  {
432  lambdaPf[pFacei] =
433  min(lambdaPf[pFacei], lambdap[pfCelli]);
434  }
435  else
436  {
437  lambdaPf[pFacei] =
438  min(lambdaPf[pFacei], lambdam[pfCelli]);
439  }
440  }
441  }
442  else
443  {
444  const labelList& pFaceCells =
445  mesh.boundary()[patchi].faceCells();
446  const scalarField& phiPf = phiBf[patchi];
447 
448  forAll(lambdaPf, pFacei)
449  {
450  // Limit outlet faces only
451  if ((phiPf[pFacei] + phiCorrfPf[pFacei]) > SMALL*SMALL)
452  {
453  const label pfCelli = pFaceCells[pFacei];
454 
455  if (phiCorrfPf[pFacei] > 0)
456  {
457  lambdaPf[pFacei] =
458  min(lambdaPf[pFacei], lambdap[pfCelli]);
459  }
460  else
461  {
462  lambdaPf[pFacei] =
463  min(lambdaPf[pFacei], lambdam[pfCelli]);
464  }
465  }
466  }
467  }
468  }
469 
470  syncTools::syncFaceList(mesh, allLambda, minEqOp<scalar>());
471  }
472 }
473 
474 
475 template<class RdeltaTType, class RhoType, class SpType, class SuType>
477 (
478  const RdeltaTType& rDeltaT,
479  const RhoType& rho,
480  const volScalarField& psi,
481  const surfaceScalarField& phi,
482  surfaceScalarField& phiCorr,
483  const SpType& Sp,
484  const SuType& Su,
485  const scalar psiMax,
486  const scalar psiMin
487 )
488 {
489  const fvMesh& mesh = psi.mesh();
490 
491  scalarField allLambda(mesh.nFaces(), 1.0);
492 
494  (
495  IOobject
496  (
497  "lambda",
498  mesh.time().timeName(),
499  mesh,
500  IOobject::NO_READ,
501  IOobject::NO_WRITE,
502  false
503  ),
504  mesh,
505  dimless,
506  allLambda,
507  false // Use slices for the couples
508  );
509 
511  (
512  allLambda,
513  rDeltaT,
514  rho,
515  psi,
516  phi,
517  phiCorr,
518  Sp,
519  Su,
520  psiMax,
521  psiMin
522  );
523 
524  phiCorr *= lambda;
525 }
526 
527 
528 // ************************************************************************* //
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const word & name() const
Return name.
Definition: IOobject.H:291
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
bool moving() const
Is mesh moving.
Definition: polyMesh.H:496
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:315
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
label nFaces() const
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:342
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:288
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
dynamicFvMesh & mesh
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:432
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:328
label readLabel(Istream &is)
Definition: label.H:64
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
const Mesh & mesh() const
Return mesh.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-5.0/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
void limitCorr(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
const volScalarField & psi
void limiterCorr(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
A class for managing temporary objects.
Definition: PtrList.H:53
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545