MULESTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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>
83 (
84  const RhoType& rho,
85  volScalarField& psi,
86  const surfaceScalarField& phiPsi
87 )
88 {
89  explicitSolve(rho, psi, phiPsi, zeroField(), zeroField());
90 }
91 
92 
93 template<class RhoType, class SpType, class SuType>
95 (
96  const RhoType& rho,
97  volScalarField& psi,
98  const surfaceScalarField& phiPsi,
99  const SpType& Sp,
100  const SuType& Su
101 )
102 {
103  const fvMesh& mesh = psi.mesh();
104 
105  if (fv::localEulerDdt::enabled(mesh))
106  {
107  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
108  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
109  }
110  else
111  {
112  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
113  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
114  }
115 }
116 
117 
118 template<class RhoType, class PsiMaxType, class PsiMinType>
120 (
121  const RhoType& rho,
122  volScalarField& psi,
123  const surfaceScalarField& phiBD,
124  surfaceScalarField& phiPsi,
125  const PsiMaxType& psiMax,
126  const PsiMinType& psiMin
127 )
128 {
130  (
131  rho,
132  psi,
133  phiBD,
134  phiPsi,
135  zeroField(),
136  zeroField(),
137  psiMax,
138  psiMin
139  );
140 }
141 
142 
143 template
144 <
145  class RhoType,
146  class SpType,
147  class SuType,
148  class PsiMaxType,
149  class PsiMinType
150 >
152 (
153  const RhoType& rho,
154  volScalarField& psi,
155  const surfaceScalarField& phi,
156  surfaceScalarField& phiPsi,
157  const SpType& Sp,
158  const SuType& Su,
159  const PsiMaxType& psiMax,
160  const PsiMinType& psiMin
161 )
162 {
163  const fvMesh& mesh = psi.mesh();
164 
166 
167  if (fv::localEulerDdt::enabled(mesh))
168  {
169  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
170  limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, false);
171  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
172  }
173  else
174  {
175  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
176  limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, false);
177  explicitSolve(rDeltaT, rho, psi, phiPsi, Sp, Su);
178  }
179 }
180 
181 
182 template
183 <
184  class RdeltaTType,
185  class RhoType,
186  class SpType,
187  class SuType,
188  class PsiMaxType,
189  class PsiMinType
190 >
192 (
193  scalarField& allLambda,
194  const RdeltaTType& rDeltaT,
195  const RhoType& rho,
196  const volScalarField& psi,
197  const surfaceScalarField& phiBD,
198  const surfaceScalarField& phiCorr,
199  const SpType& Sp,
200  const SuType& Su,
201  const PsiMaxType& psiMax,
202  const PsiMinType& psiMin
203 )
204 {
205  const scalarField& psiIf = psi;
206  const volScalarField::Boundary& psiBf = psi.boundaryField();
207 
208  const fvMesh& mesh = psi.mesh();
209 
210  const dictionary& MULEScontrols = mesh.solverDict(psi.name());
211 
212  const label nLimiterIter
213  (
214  MULEScontrols.lookupOrDefault<label>("nLimiterIter", 3)
215  );
216 
217  const scalar smoothLimiter
218  (
219  MULEScontrols.lookupOrDefault<scalar>("smoothLimiter", 0)
220  );
221 
222  const scalar extremaCoeff
223  (
224  MULEScontrols.lookupOrDefault<scalar>("extremaCoeff", 0)
225  );
226 
227  const scalarField& psi0 = psi.oldTime();
228 
229  const labelUList& owner = mesh.owner();
230  const labelUList& neighb = mesh.neighbour();
231  tmp<volScalarField::Internal> tVsc = mesh.Vsc();
232  const scalarField& V = tVsc();
233 
234  const scalarField& phiBDIf = phiBD;
235  const surfaceScalarField::Boundary& phiBDBf =
236  phiBD.boundaryField();
237 
238  const scalarField& phiCorrIf = phiCorr;
239  const surfaceScalarField::Boundary& phiCorrBf =
240  phiCorr.boundaryField();
241 
243  (
244  IOobject
245  (
246  "lambda",
247  mesh.time().timeName(),
248  mesh,
249  IOobject::NO_READ,
250  IOobject::NO_WRITE,
251  false
252  ),
253  mesh,
254  dimless,
255  allLambda,
256  false // Use slices for the couples
257  );
258 
259  scalarField& lambdaIf = lambda;
260  surfaceScalarField::Boundary& lambdaBf = lambda.boundaryFieldRef();
261 
262  scalarField psiMaxn(psiIf.size());
263  scalarField psiMinn(psiIf.size());
264 
265  psiMaxn = psiMin;
266  psiMinn = psiMax;
267 
268  scalarField sumPhiBD(psiIf.size(), 0.0);
269 
270  scalarField sumPhip(psiIf.size(), 0.0);
271  scalarField mSumPhim(psiIf.size(), 0.0);
272 
273  forAll(phiCorrIf, facei)
274  {
275  const label own = owner[facei];
276  const label nei = neighb[facei];
277 
278  psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
279  psiMinn[own] = min(psiMinn[own], psiIf[nei]);
280 
281  psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
282  psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
283 
284  sumPhiBD[own] += phiBDIf[facei];
285  sumPhiBD[nei] -= phiBDIf[facei];
286 
287  const scalar phiCorrf = phiCorrIf[facei];
288 
289  if (phiCorrf > 0)
290  {
291  sumPhip[own] += phiCorrf;
292  mSumPhim[nei] += phiCorrf;
293  }
294  else
295  {
296  mSumPhim[own] -= phiCorrf;
297  sumPhip[nei] -= phiCorrf;
298  }
299  }
300 
301  forAll(phiCorrBf, patchi)
302  {
303  const fvPatchScalarField& psiPf = psiBf[patchi];
304  const scalarField& phiBDPf = phiBDBf[patchi];
305  const scalarField& phiCorrPf = phiCorrBf[patchi];
306 
307  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
308 
309  if (psiPf.coupled())
310  {
311  const scalarField psiPNf(psiPf.patchNeighbourField());
312 
313  forAll(phiCorrPf, pFacei)
314  {
315  const label pfCelli = pFaceCells[pFacei];
316 
317  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
318  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
319  }
320  }
321  else if (psiPf.fixesValue())
322  {
323  forAll(phiCorrPf, pFacei)
324  {
325  const label pfCelli = pFaceCells[pFacei];
326 
327  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
328  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
329  }
330  }
331  else
332  {
333  forAll(phiCorrPf, pFacei)
334  {
335  const label pfCelli = pFaceCells[pFacei];
336 
337  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiMax[pfCelli]);
338  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiMin[pfCelli]);
339  }
340  }
341 
342  forAll(phiCorrPf, pFacei)
343  {
344  const label pfCelli = pFaceCells[pFacei];
345 
346  sumPhiBD[pfCelli] += phiBDPf[pFacei];
347 
348  const scalar phiCorrf = phiCorrPf[pFacei];
349 
350  if (phiCorrf > 0)
351  {
352  sumPhip[pfCelli] += phiCorrf;
353  }
354  else
355  {
356  mSumPhim[pfCelli] -= phiCorrf;
357  }
358  }
359  }
360 
361  psiMaxn = min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
362  psiMinn = max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
363 
364  if (smoothLimiter > small)
365  {
366  psiMaxn =
367  min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
368  psiMinn =
369  max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
370  }
371 
372  if (mesh.moving())
373  {
375 
376  psiMaxn =
377  V
378  *(
379  (rho.field()*rDeltaT - Sp.field())*psiMaxn
380  - Su.field()
381  )
382  - (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
383  + sumPhiBD;
384 
385  psiMinn =
386  V
387  *(
388  Su.field()
389  - (rho.field()*rDeltaT - Sp.field())*psiMinn
390  )
391  + (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
392  - sumPhiBD;
393  }
394  else
395  {
396  psiMaxn =
397  V
398  *(
399  (rho.field()*rDeltaT - Sp.field())*psiMaxn
400  - Su.field()
401  - (rho.oldTime().field()*rDeltaT)*psi0
402  )
403  + sumPhiBD;
404 
405  psiMinn =
406  V
407  *(
408  Su.field()
409  - (rho.field()*rDeltaT - Sp.field())*psiMinn
410  + (rho.oldTime().field()*rDeltaT)*psi0
411  )
412  - sumPhiBD;
413  }
414 
415  scalarField sumlPhip(psiIf.size());
416  scalarField mSumlPhim(psiIf.size());
417 
418  for (int j=0; j<nLimiterIter; j++)
419  {
420  sumlPhip = 0;
421  mSumlPhim = 0;
422 
423  forAll(lambdaIf, facei)
424  {
425  const label own = owner[facei];
426  const label nei = neighb[facei];
427 
428  scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
429 
430  if (lambdaPhiCorrf > 0)
431  {
432  sumlPhip[own] += lambdaPhiCorrf;
433  mSumlPhim[nei] += lambdaPhiCorrf;
434  }
435  else
436  {
437  mSumlPhim[own] -= lambdaPhiCorrf;
438  sumlPhip[nei] -= lambdaPhiCorrf;
439  }
440  }
441 
442  forAll(lambdaBf, patchi)
443  {
444  scalarField& lambdaPf = lambdaBf[patchi];
445  const scalarField& phiCorrfPf = phiCorrBf[patchi];
446 
447  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
448 
449  forAll(lambdaPf, pFacei)
450  {
451  const label pfCelli = pFaceCells[pFacei];
452  const scalar lambdaPhiCorrf =
453  lambdaPf[pFacei]*phiCorrfPf[pFacei];
454 
455  if (lambdaPhiCorrf > 0)
456  {
457  sumlPhip[pfCelli] += lambdaPhiCorrf;
458  }
459  else
460  {
461  mSumlPhim[pfCelli] -= lambdaPhiCorrf;
462  }
463  }
464  }
465 
466  forAll(sumlPhip, celli)
467  {
468  sumlPhip[celli] =
469  max(min
470  (
471  (sumlPhip[celli] + psiMaxn[celli])
472  /(mSumPhim[celli] + rootVSmall),
473  1.0), 0.0
474  );
475 
476  mSumlPhim[celli] =
477  max(min
478  (
479  (mSumlPhim[celli] + psiMinn[celli])
480  /(sumPhip[celli] + rootVSmall),
481  1.0), 0.0
482  );
483  }
484 
485  const scalarField& lambdam = sumlPhip;
486  const scalarField& lambdap = mSumlPhim;
487 
488  forAll(lambdaIf, facei)
489  {
490  if (phiCorrIf[facei] > 0)
491  {
492  lambdaIf[facei] = min
493  (
494  lambdaIf[facei],
495  min(lambdap[owner[facei]], lambdam[neighb[facei]])
496  );
497  }
498  else
499  {
500  lambdaIf[facei] = min
501  (
502  lambdaIf[facei],
503  min(lambdam[owner[facei]], lambdap[neighb[facei]])
504  );
505  }
506  }
507 
508  forAll(lambdaBf, patchi)
509  {
510  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
511  const scalarField& phiCorrfPf = phiCorrBf[patchi];
512  const fvPatchScalarField& psiPf = psiBf[patchi];
513 
514  if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
515  {
516  lambdaPf = 0;
517  }
518  else if (psiPf.coupled())
519  {
520  const labelList& pFaceCells =
521  mesh.boundary()[patchi].faceCells();
522 
523  forAll(lambdaPf, pFacei)
524  {
525  const label pfCelli = pFaceCells[pFacei];
526 
527  if (phiCorrfPf[pFacei] > 0)
528  {
529  lambdaPf[pFacei] =
530  min(lambdaPf[pFacei], lambdap[pfCelli]);
531  }
532  else
533  {
534  lambdaPf[pFacei] =
535  min(lambdaPf[pFacei], lambdam[pfCelli]);
536  }
537  }
538  }
539  }
540 
541  syncTools::syncFaceList(mesh, allLambda, minEqOp<scalar>());
542  }
543 }
544 
545 
546 template
547 <
548  class RdeltaTType,
549  class RhoType,
550  class SpType,
551  class SuType,
552  class PsiMaxType,
553  class PsiMinType
554 >
556 (
557  const RdeltaTType& rDeltaT,
558  const RhoType& rho,
559  const volScalarField& psi,
560  const surfaceScalarField& phi,
561  surfaceScalarField& phiPsi,
562  const SpType& Sp,
563  const SuType& Su,
564  const PsiMaxType& psiMax,
565  const PsiMinType& psiMin,
566  const bool returnCorr
567 )
568 {
569  const fvMesh& mesh = psi.mesh();
570 
571  surfaceScalarField phiBD(upwind<scalar>(psi.mesh(), phi).flux(psi));
572 
573  surfaceScalarField::Boundary& phiBDBf = phiBD.boundaryFieldRef();
574  const surfaceScalarField::Boundary& phiPsiBf = phiPsi.boundaryField();
575 
576  forAll(phiBDBf, patchi)
577  {
578  fvsPatchScalarField& phiBDPf = phiBDBf[patchi];
579 
580  if (!phiBDPf.coupled())
581  {
582  phiBDPf = phiPsiBf[patchi];
583  }
584  }
585 
586  surfaceScalarField& phiCorr = phiPsi;
587  phiCorr -= phiBD;
588 
589  scalarField allLambda(mesh.nFaces(), 1.0);
590 
592  (
593  IOobject
594  (
595  "lambda",
596  mesh.time().timeName(),
597  mesh,
598  IOobject::NO_READ,
599  IOobject::NO_WRITE,
600  false
601  ),
602  mesh,
603  dimless,
604  allLambda,
605  false // Use slices for the couples
606  );
607 
608  limiter
609  (
610  allLambda,
611  rDeltaT,
612  rho,
613  psi,
614  phiBD,
615  phiCorr,
616  Sp,
617  Su,
618  psiMax,
619  psiMin
620  );
621 
622  if (returnCorr)
623  {
624  phiCorr *= lambda;
625  }
626  else
627  {
628  phiPsi = phiBD + lambda*phiCorr;
629  }
630 }
631 
632 
633 template
634 <
635  class RhoType,
636  class SpType,
637  class SuType,
638  class PsiMaxType,
639  class PsiMinType
640 >
642 (
643  const RhoType& rho,
644  const volScalarField& psi,
645  const surfaceScalarField& phi,
646  surfaceScalarField& phiPsi,
647  const SpType& Sp,
648  const SuType& Su,
649  const PsiMaxType& psiMax,
650  const PsiMinType& psiMin,
651  const bool rtnCorr
652 )
653 {
654  const fvMesh& mesh = psi.mesh();
655 
656  if (fv::localEulerDdt::enabled(mesh))
657  {
658  const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
659  limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
660  }
661  else
662  {
663  const scalar rDeltaT = 1.0/mesh.time().deltaTValue();
664  limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
665  }
666 }
667 
668 
669 template<class SurfaceScalarFieldList>
670 void Foam::MULES::limitSum(SurfaceScalarFieldList& phiPsiCorrs)
671 {
672  {
673  UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
674  forAll(phiPsiCorrs, phasei)
675  {
676  phiPsiCorrsInternal.set(phasei, &phiPsiCorrs[phasei]);
677  }
678 
679  limitSum(phiPsiCorrsInternal);
680  }
681 
682  const surfaceScalarField::Boundary& bfld =
683  phiPsiCorrs[0].boundaryField();
684 
685  forAll(bfld, patchi)
686  {
687  if (bfld[patchi].coupled())
688  {
689  UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
690  forAll(phiPsiCorrs, phasei)
691  {
692  phiPsiCorrsPatch.set
693  (
694  phasei,
695  &phiPsiCorrs[phasei].boundaryFieldRef()[patchi]
696  );
697  }
698 
699  limitSum(phiPsiCorrsPatch);
700  }
701  }
702 }
703 
704 
705 template<class SurfaceScalarFieldList>
707 (
708  const SurfaceScalarFieldList& alphas,
709  SurfaceScalarFieldList& phiPsiCorrs,
710  const labelHashSet& fixed
711 )
712 {
713  {
714  UPtrList<const scalarField> alphasInternal(alphas.size());
715  forAll(alphas, phasei)
716  {
717  alphasInternal.set(phasei, &alphas[phasei]);
718  }
719  UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
720  forAll(phiPsiCorrs, phasei)
721  {
722  phiPsiCorrsInternal.set(phasei, &phiPsiCorrs[phasei]);
723  }
724 
725  limitSum(alphasInternal, phiPsiCorrsInternal, fixed);
726  }
727 
728  const surfaceScalarField::Boundary& bfld =
729  phiPsiCorrs[0].boundaryField();
730 
731  forAll(bfld, patchi)
732  {
733  if (bfld[patchi].coupled())
734  {
735  UPtrList<const scalarField> alphasPatch(alphas.size());
736  forAll(alphas, phasei)
737  {
738  alphasPatch.set
739  (
740  phasei,
741  &alphas[phasei].boundaryField()[patchi]
742  );
743  }
744  UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
745  forAll(phiPsiCorrs, phasei)
746  {
747  phiPsiCorrsPatch.set
748  (
749  phasei,
750  &phiPsiCorrs[phasei].boundaryFieldRef()[patchi]
751  );
752  }
753 
754  limitSum(alphasPatch, phiPsiCorrsPatch, fixed);
755  }
756  }
757 }
758 
759 
760 // ************************************************************************* //
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
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
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
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:297
surfaceScalarField & phi
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
bool moving() const
Is mesh moving.
Definition: polyMesh.H:502
label phasei
Definition: pEqn.H:27
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:256
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
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
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:30
virtual bool coupled() const
Return true if this patch field is coupled.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:50
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 PsiMaxType &psiMax, const PsiMinType &psiMin)
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
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
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
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycl old-time cell volumes.
pressureControl limit(p)
const Mesh & mesh() const
Return mesh.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
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.
Upwind differencing scheme class.
Definition: upwind.H:51
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
MULES: Multidimensional universal limiter for explicit solution.
const volScalarField & psi
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
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
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545