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-2016 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  label nLimiterIter
155  (
156  readLabel(MULEScontrols.lookup("nLimiterIter"))
157  );
158 
159  scalar smoothLimiter
160  (
161  MULEScontrols.lookupOrDefault<scalar>("smoothLimiter", 0)
162  );
163 
164  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(), VSMALL);
206  scalarField mSumPhim(psiIf.size(), VSMALL);
207 
208  forAll(phiCorrIf, facei)
209  {
210  label own = owner[facei];
211  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  scalar phiCorrf = phiCorrIf[facei];
220 
221  if (phiCorrf > 0.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
253  {
254  forAll(phiCorrPf, pFacei)
255  {
256  label pfCelli = pFaceCells[pFacei];
257 
258  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
259  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
260  }
261  }
262 
263  forAll(phiCorrPf, pFacei)
264  {
265  label pfCelli = pFaceCells[pFacei];
266 
267  scalar phiCorrf = phiCorrPf[pFacei];
268 
269  if (phiCorrf > 0.0)
270  {
271  sumPhip[pfCelli] += phiCorrf;
272  }
273  else
274  {
275  mSumPhim[pfCelli] -= phiCorrf;
276  }
277  }
278  }
279 
280  psiMaxn = min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
281  psiMinn = max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
282 
283  if (smoothLimiter > SMALL)
284  {
285  psiMaxn =
286  min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
287  psiMinn =
288  max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
289  }
290 
291  psiMaxn =
292  V
293  *(
294  (rho.field()*rDeltaT - Sp.field())*psiMaxn
295  - Su.field()
296  - rho.field()*psi.primitiveField()*rDeltaT
297  );
298 
299  psiMinn =
300  V
301  *(
302  Su.field()
303  - (rho.field()*rDeltaT - Sp.field())*psiMinn
304  + rho.field()*psi.primitiveField()*rDeltaT
305  );
306 
307  scalarField sumlPhip(psiIf.size());
308  scalarField mSumlPhim(psiIf.size());
309 
310  for (int j=0; j<nLimiterIter; j++)
311  {
312  sumlPhip = 0.0;
313  mSumlPhim = 0.0;
314 
315  forAll(lambdaIf, facei)
316  {
317  label own = owner[facei];
318  label nei = neighb[facei];
319 
320  scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
321 
322  if (lambdaPhiCorrf > 0.0)
323  {
324  sumlPhip[own] += lambdaPhiCorrf;
325  mSumlPhim[nei] += lambdaPhiCorrf;
326  }
327  else
328  {
329  mSumlPhim[own] -= lambdaPhiCorrf;
330  sumlPhip[nei] -= lambdaPhiCorrf;
331  }
332  }
333 
334  forAll(lambdaBf, patchi)
335  {
336  scalarField& lambdaPf = lambdaBf[patchi];
337  const scalarField& phiCorrfPf = phiCorrBf[patchi];
338 
339  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
340 
341  forAll(lambdaPf, pFacei)
342  {
343  label pfCelli = pFaceCells[pFacei];
344 
345  scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
346 
347  if (lambdaPhiCorrf > 0.0)
348  {
349  sumlPhip[pfCelli] += lambdaPhiCorrf;
350  }
351  else
352  {
353  mSumlPhim[pfCelli] -= lambdaPhiCorrf;
354  }
355  }
356  }
357 
358  forAll(sumlPhip, celli)
359  {
360  sumlPhip[celli] =
361  max(min
362  (
363  (sumlPhip[celli] + psiMaxn[celli])
364  /(mSumPhim[celli] - SMALL),
365  1.0), 0.0
366  );
367 
368  mSumlPhim[celli] =
369  max(min
370  (
371  (mSumlPhim[celli] + psiMinn[celli])
372  /(sumPhip[celli] + SMALL),
373  1.0), 0.0
374  );
375  }
376 
377  const scalarField& lambdam = sumlPhip;
378  const scalarField& lambdap = mSumlPhim;
379 
380  forAll(lambdaIf, facei)
381  {
382  if (phiCorrIf[facei] > 0.0)
383  {
384  lambdaIf[facei] = min
385  (
386  lambdaIf[facei],
387  min(lambdap[owner[facei]], lambdam[neighb[facei]])
388  );
389  }
390  else
391  {
392  lambdaIf[facei] = min
393  (
394  lambdaIf[facei],
395  min(lambdam[owner[facei]], lambdap[neighb[facei]])
396  );
397  }
398  }
399 
400 
401  forAll(lambdaBf, patchi)
402  {
403  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
404  const scalarField& phiCorrfPf = phiCorrBf[patchi];
405  const fvPatchScalarField& psiPf = psiBf[patchi];
406 
407  if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
408  {
409  lambdaPf = 0;
410  }
411  else if (psiPf.coupled())
412  {
413  const labelList& pFaceCells =
414  mesh.boundary()[patchi].faceCells();
415 
416  forAll(lambdaPf, pFacei)
417  {
418  label pfCelli = pFaceCells[pFacei];
419 
420  if (phiCorrfPf[pFacei] > 0.0)
421  {
422  lambdaPf[pFacei] =
423  min(lambdaPf[pFacei], lambdap[pfCelli]);
424  }
425  else
426  {
427  lambdaPf[pFacei] =
428  min(lambdaPf[pFacei], lambdam[pfCelli]);
429  }
430  }
431  }
432  else
433  {
434  const labelList& pFaceCells =
435  mesh.boundary()[patchi].faceCells();
436  const scalarField& phiPf = phiBf[patchi];
437 
438  forAll(lambdaPf, pFacei)
439  {
440  // Limit outlet faces only
441  if (phiPf[pFacei] > SMALL*SMALL)
442  {
443  label pfCelli = pFaceCells[pFacei];
444 
445  if (phiCorrfPf[pFacei] > 0.0)
446  {
447  lambdaPf[pFacei] =
448  min(lambdaPf[pFacei], lambdap[pfCelli]);
449  }
450  else
451  {
452  lambdaPf[pFacei] =
453  min(lambdaPf[pFacei], lambdam[pfCelli]);
454  }
455  }
456  }
457  }
458  }
459 
460  syncTools::syncFaceList(mesh, allLambda, minEqOp<scalar>());
461  }
462 }
463 
464 
465 template<class RdeltaTType, class RhoType, class SpType, class SuType>
467 (
468  const RdeltaTType& rDeltaT,
469  const RhoType& rho,
470  const volScalarField& psi,
471  const surfaceScalarField& phi,
472  surfaceScalarField& phiCorr,
473  const SpType& Sp,
474  const SuType& Su,
475  const scalar psiMax,
476  const scalar psiMin
477 )
478 {
479  const fvMesh& mesh = psi.mesh();
480 
481  scalarField allLambda(mesh.nFaces(), 1.0);
482 
484  (
485  IOobject
486  (
487  "lambda",
488  mesh.time().timeName(),
489  mesh,
490  IOobject::NO_READ,
491  IOobject::NO_WRITE,
492  false
493  ),
494  mesh,
495  dimless,
496  allLambda,
497  false // Use slices for the couples
498  );
499 
501  (
502  allLambda,
503  rDeltaT,
504  rho,
505  psi,
506  phi,
507  phiCorr,
508  Sp,
509  Su,
510  psiMax,
511  psiMin
512  );
513 
514  phiCorr *= lambda;
515 }
516 
517 
518 // ************************************************************************* //
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
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:327
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
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 > &)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
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:65
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:715
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:431
dynamicFvMesh & mesh
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:288
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
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:60
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
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-4.1/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
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label nFaces() const
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
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const Mesh & mesh() const
Return mesh.
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
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:342
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
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:54
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:91
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
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.
const word & name() const
Return name.
Definition: IOobject.H:260
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451