MULESTemplates.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) 2011-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 "MULES.H"
27 #include "upwind.H"
28 #include "fvcSurfaceIntegrate.H"
29 #include "localEulerDdtScheme.H"
30 #include "slicedSurfaceFields.H"
31 #include "wedgeFvPatch.H"
32 #include "syncTools.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 template<class RdeltaTType, class RhoType, class SpType, class SuType>
38 (
39  const RdeltaTType& rDeltaT,
40  const RhoType& rho,
41  volScalarField& psi,
42  const surfaceScalarField& phiPsi,
43  const SpType& Sp,
44  const SuType& Su
45 )
46 {
47  Info<< "MULES: Solving for " << psi.name() << endl;
48 
49  const fvMesh& mesh = psi.mesh();
50 
51  scalarField& psiIf = psi;
52  const scalarField& psi0 = psi.oldTime();
53 
54  psiIf = 0.0;
55  fvc::surfaceIntegrate(psiIf, phiPsi);
56 
57  if (mesh.moving())
58  {
59  psiIf =
60  (
61  mesh.Vsc0()().field()*rho.oldTime().field()
62  *psi0*rDeltaT/mesh.Vsc()().field()
63  + Su.field()
64  - psiIf
65  )/(rho.field()*rDeltaT - Sp.field());
66  }
67  else
68  {
69  psiIf =
70  (
71  rho.oldTime().field()*psi0*rDeltaT
72  + Su.field()
73  - psiIf
74  )/(rho.field()*rDeltaT - Sp.field());
75  }
76 
78 }
79 
80 
81 template<class RhoType, class SpType, class SuType>
83 (
84  const RhoType& rho,
85  volScalarField& psi,
86  const surfaceScalarField& phiPsi,
87  const SpType& Sp,
88  const SuType& Su
89 )
90 {
91  const fvMesh& mesh = psi.mesh();
92 
93  if (fv::localEulerDdt::enabled(mesh))
94  {
95  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
96  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
97  }
98  else
99  {
100  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
101  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
102  }
103 }
104 
105 
106 template<class RhoType, class SpType, class SuType>
108 (
109  const RhoType& rho,
110  volScalarField& psi,
111  const surfaceScalarField& phi,
112  surfaceScalarField& phiPsi,
113  const SpType& Sp,
114  const SuType& Su,
115  const scalar psiMax,
116  const scalar psiMin
117 )
118 {
119  const fvMesh& mesh = psi.mesh();
120 
122 
123  if (fv::localEulerDdt::enabled(mesh))
124  {
125  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
126 
127  limit
128  (
129  rDeltaT,
130  rho,
131  psi,
132  phi,
133  phiPsi,
134  Sp,
135  Su,
136  psiMax,
137  psiMin,
138  false
139  );
140 
141  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
142  }
143  else
144  {
145  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
146 
147  limit
148  (
149  rDeltaT,
150  rho,
151  psi,
152  phi,
153  phiPsi,
154  Sp,
155  Su,
156  psiMax,
157  psiMin,
158  false
159  );
160 
161  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
162  }
163 }
164 
165 
166 template<class RdeltaTType, class RhoType, class SpType, class SuType>
168 (
169  scalarField& allLambda,
170  const RdeltaTType& rDeltaT,
171  const RhoType& rho,
172  const volScalarField& psi,
173  const surfaceScalarField& phiBD,
174  const surfaceScalarField& phiCorr,
175  const SpType& Sp,
176  const SuType& Su,
177  const scalar psiMax,
178  const scalar psiMin
179 )
180 {
181  const scalarField& psiIf = psi;
182  const volScalarField::GeometricBoundaryField& psiBf = psi.boundaryField();
183 
184  const fvMesh& mesh = psi.mesh();
185 
186  const dictionary& MULEScontrols = mesh.solverDict(psi.name());
187 
188  label nLimiterIter
189  (
190  MULEScontrols.lookupOrDefault<label>("nLimiterIter", 3)
191  );
192 
193  scalar smoothLimiter
194  (
195  MULEScontrols.lookupOrDefault<scalar>("smoothLimiter", 0)
196  );
197 
198  const scalarField& psi0 = psi.oldTime();
199 
200  const labelUList& owner = mesh.owner();
201  const labelUList& neighb = mesh.neighbour();
203  const scalarField& V = tVsc();
204 
205  const scalarField& phiBDIf = phiBD;
206  const surfaceScalarField::GeometricBoundaryField& phiBDBf =
207  phiBD.boundaryField();
208 
209  const scalarField& phiCorrIf = phiCorr;
210  const surfaceScalarField::GeometricBoundaryField& phiCorrBf =
211  phiCorr.boundaryField();
212 
214  (
215  IOobject
216  (
217  "lambda",
218  mesh.time().timeName(),
219  mesh,
220  IOobject::NO_READ,
221  IOobject::NO_WRITE,
222  false
223  ),
224  mesh,
225  dimless,
226  allLambda,
227  false // Use slices for the couples
228  );
229 
230  scalarField& lambdaIf = lambda;
231  surfaceScalarField::GeometricBoundaryField& lambdaBf =
232  lambda.boundaryField();
233 
234  scalarField psiMaxn(psiIf.size(), psiMin);
235  scalarField psiMinn(psiIf.size(), psiMax);
236 
237  scalarField sumPhiBD(psiIf.size(), 0.0);
238 
239  scalarField sumPhip(psiIf.size(), VSMALL);
240  scalarField mSumPhim(psiIf.size(), VSMALL);
241 
242  forAll(phiCorrIf, facei)
243  {
244  label own = owner[facei];
245  label nei = neighb[facei];
246 
247  psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
248  psiMinn[own] = min(psiMinn[own], psiIf[nei]);
249 
250  psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
251  psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
252 
253  sumPhiBD[own] += phiBDIf[facei];
254  sumPhiBD[nei] -= phiBDIf[facei];
255 
256  scalar phiCorrf = phiCorrIf[facei];
257 
258  if (phiCorrf > 0.0)
259  {
260  sumPhip[own] += phiCorrf;
261  mSumPhim[nei] += phiCorrf;
262  }
263  else
264  {
265  mSumPhim[own] -= phiCorrf;
266  sumPhip[nei] -= phiCorrf;
267  }
268  }
269 
270  forAll(phiCorrBf, patchi)
271  {
272  const fvPatchScalarField& psiPf = psiBf[patchi];
273  const scalarField& phiBDPf = phiBDBf[patchi];
274  const scalarField& phiCorrPf = phiCorrBf[patchi];
275 
276  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
277 
278  if (psiPf.coupled())
279  {
280  const scalarField psiPNf(psiPf.patchNeighbourField());
281 
282  forAll(phiCorrPf, pFacei)
283  {
284  label pfCelli = pFaceCells[pFacei];
285 
286  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
287  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
288  }
289  }
290  else
291  {
292  forAll(phiCorrPf, pFacei)
293  {
294  label pfCelli = pFaceCells[pFacei];
295 
296  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
297  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
298  }
299  }
300 
301  forAll(phiCorrPf, pFacei)
302  {
303  label pfCelli = pFaceCells[pFacei];
304 
305  sumPhiBD[pfCelli] += phiBDPf[pFacei];
306 
307  scalar phiCorrf = phiCorrPf[pFacei];
308 
309  if (phiCorrf > 0.0)
310  {
311  sumPhip[pfCelli] += phiCorrf;
312  }
313  else
314  {
315  mSumPhim[pfCelli] -= phiCorrf;
316  }
317  }
318  }
319 
320  psiMaxn = min(psiMaxn, psiMax);
321  psiMinn = max(psiMinn, psiMin);
322 
323  if (smoothLimiter > SMALL)
324  {
325  psiMaxn =
326  min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
327  psiMinn =
328  max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
329  }
330 
331  if (mesh.moving())
332  {
334 
335  psiMaxn =
336  V
337  *(
338  (rho.field()*rDeltaT - Sp.field())*psiMaxn
339  - Su.field()
340  )
341  - (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
342  + sumPhiBD;
343 
344  psiMinn =
345  V
346  *(
347  Su.field()
348  - (rho.field()*rDeltaT - Sp.field())*psiMinn
349  )
350  + (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
351  - sumPhiBD;
352  }
353  else
354  {
355  psiMaxn =
356  V
357  *(
358  (rho.field()*rDeltaT - Sp.field())*psiMaxn
359  - Su.field()
360  - (rho.oldTime().field()*rDeltaT)*psi0
361  )
362  + sumPhiBD;
363 
364  psiMinn =
365  V
366  *(
367  Su.field()
368  - (rho.field()*rDeltaT - Sp.field())*psiMinn
369  + (rho.oldTime().field()*rDeltaT)*psi0
370  )
371  - sumPhiBD;
372  }
373 
374  scalarField sumlPhip(psiIf.size());
375  scalarField mSumlPhim(psiIf.size());
376 
377  for (int j=0; j<nLimiterIter; j++)
378  {
379  sumlPhip = 0.0;
380  mSumlPhim = 0.0;
381 
382  forAll(lambdaIf, facei)
383  {
384  label own = owner[facei];
385  label nei = neighb[facei];
386 
387  scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
388 
389  if (lambdaPhiCorrf > 0.0)
390  {
391  sumlPhip[own] += lambdaPhiCorrf;
392  mSumlPhim[nei] += lambdaPhiCorrf;
393  }
394  else
395  {
396  mSumlPhim[own] -= lambdaPhiCorrf;
397  sumlPhip[nei] -= lambdaPhiCorrf;
398  }
399  }
400 
401  forAll(lambdaBf, patchi)
402  {
403  scalarField& lambdaPf = lambdaBf[patchi];
404  const scalarField& phiCorrfPf = phiCorrBf[patchi];
405 
406  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
407 
408  forAll(lambdaPf, pFacei)
409  {
410  label pfCelli = pFaceCells[pFacei];
411 
412  scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
413 
414  if (lambdaPhiCorrf > 0.0)
415  {
416  sumlPhip[pfCelli] += lambdaPhiCorrf;
417  }
418  else
419  {
420  mSumlPhim[pfCelli] -= lambdaPhiCorrf;
421  }
422  }
423  }
424 
425  forAll(sumlPhip, celli)
426  {
427  sumlPhip[celli] =
428  max(min
429  (
430  (sumlPhip[celli] + psiMaxn[celli])
431  /(mSumPhim[celli] - SMALL),
432  1.0), 0.0
433  );
434 
435  mSumlPhim[celli] =
436  max(min
437  (
438  (mSumlPhim[celli] + psiMinn[celli])
439  /(sumPhip[celli] + SMALL),
440  1.0), 0.0
441  );
442  }
443 
444  const scalarField& lambdam = sumlPhip;
445  const scalarField& lambdap = mSumlPhim;
446 
447  forAll(lambdaIf, facei)
448  {
449  if (phiCorrIf[facei] > 0.0)
450  {
451  lambdaIf[facei] = min
452  (
453  lambdaIf[facei],
454  min(lambdap[owner[facei]], lambdam[neighb[facei]])
455  );
456  }
457  else
458  {
459  lambdaIf[facei] = min
460  (
461  lambdaIf[facei],
462  min(lambdam[owner[facei]], lambdap[neighb[facei]])
463  );
464  }
465  }
466 
467 
468  forAll(lambdaBf, patchi)
469  {
470  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
471  const scalarField& phiCorrfPf = phiCorrBf[patchi];
472  const fvPatchScalarField& psiPf = psiBf[patchi];
473 
474  if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
475  {
476  lambdaPf = 0;
477  }
478  else if (psiPf.coupled())
479  {
480  const labelList& pFaceCells =
481  mesh.boundary()[patchi].faceCells();
482 
483  forAll(lambdaPf, pFacei)
484  {
485  label pfCelli = pFaceCells[pFacei];
486 
487  if (phiCorrfPf[pFacei] > 0.0)
488  {
489  lambdaPf[pFacei] =
490  min(lambdaPf[pFacei], lambdap[pfCelli]);
491  }
492  else
493  {
494  lambdaPf[pFacei] =
495  min(lambdaPf[pFacei], lambdam[pfCelli]);
496  }
497  }
498  }
499  else
500  {
501  const labelList& pFaceCells =
502  mesh.boundary()[patchi].faceCells();
503  const scalarField& phiBDPf = phiBDBf[patchi];
504  const scalarField& phiCorrPf = phiCorrBf[patchi];
505 
506  forAll(lambdaPf, pFacei)
507  {
508  // Limit outlet faces only
509  if ((phiBDPf[pFacei] + phiCorrPf[pFacei]) > SMALL*SMALL)
510  {
511  label pfCelli = pFaceCells[pFacei];
512 
513  if (phiCorrfPf[pFacei] > 0.0)
514  {
515  lambdaPf[pFacei] =
516  min(lambdaPf[pFacei], lambdap[pfCelli]);
517  }
518  else
519  {
520  lambdaPf[pFacei] =
521  min(lambdaPf[pFacei], lambdam[pfCelli]);
522  }
523  }
524  }
525  }
526  }
527 
528  syncTools::syncFaceList(mesh, allLambda, minEqOp<scalar>());
529  }
530 }
531 
532 
533 template<class RdeltaTType, class RhoType, class SpType, class SuType>
535 (
536  const RdeltaTType& rDeltaT,
537  const RhoType& rho,
538  const volScalarField& psi,
539  const surfaceScalarField& phi,
540  surfaceScalarField& phiPsi,
541  const SpType& Sp,
542  const SuType& Su,
543  const scalar psiMax,
544  const scalar psiMin,
545  const bool returnCorr
546 )
547 {
548  const fvMesh& mesh = psi.mesh();
549 
550  surfaceScalarField phiBD(upwind<scalar>(psi.mesh(), phi).flux(psi));
551 
552  surfaceScalarField& phiCorr = phiPsi;
553  phiCorr -= phiBD;
554 
555  scalarField allLambda(mesh.nFaces(), 1.0);
556 
558  (
559  IOobject
560  (
561  "lambda",
562  mesh.time().timeName(),
563  mesh,
564  IOobject::NO_READ,
565  IOobject::NO_WRITE,
566  false
567  ),
568  mesh,
569  dimless,
570  allLambda,
571  false // Use slices for the couples
572  );
573 
574  limiter
575  (
576  allLambda,
577  rDeltaT,
578  rho,
579  psi,
580  phiBD,
581  phiCorr,
582  Sp,
583  Su,
584  psiMax,
585  psiMin
586  );
587 
588  if (returnCorr)
589  {
590  phiCorr *= lambda;
591  }
592  else
593  {
594  phiPsi = phiBD + lambda*phiCorr;
595  }
596 }
597 
598 
599 template<class SurfaceScalarFieldList>
600 void Foam::MULES::limitSum(SurfaceScalarFieldList& phiPsiCorrs)
601 {
602  {
603  UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
604  forAll(phiPsiCorrs, phasei)
605  {
606  phiPsiCorrsInternal.set(phasei, &phiPsiCorrs[phasei]);
607  }
608 
609  limitSum(phiPsiCorrsInternal);
610  }
611 
612  surfaceScalarField::GeometricBoundaryField& bfld =
613  phiPsiCorrs[0].boundaryField();
614 
615  forAll(bfld, patchi)
616  {
617  if (bfld[patchi].coupled())
618  {
619  UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
620  forAll(phiPsiCorrs, phasei)
621  {
622  phiPsiCorrsPatch.set
623  (
624  phasei,
625  &phiPsiCorrs[phasei].boundaryField()[patchi]
626  );
627  }
628 
629  limitSum(phiPsiCorrsPatch);
630  }
631  }
632 }
633 
634 
635 // ************************************************************************* //
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:52
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcFlux.C:45
label nFaces() const
label phasei
Definition: pEqn.H:37
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:365
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
messageStream Info
const Mesh & mesh() const
Return mesh.
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:342
MULES: Multidimensional universal limiter for explicit solution.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:411
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
scalar deltaTValue() const
Return time step value.
Definition: TimeState.H:100
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:51
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
label patchi
const volScalarField & psi
Definition: createFields.H:24
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
const word & name() const
Return name.
Definition: IOobject.H:260
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
surfaceScalarField & phi
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin, const bool returnCorr)
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
void correctBoundaryConditions()
Correct boundary field.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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
Upwind differencing scheme class.
Definition: upwind.H:51
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:288
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void limiter(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
A class for managing temporary objects.
Definition: PtrList.H:118
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycl old-time cell volumes.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:552