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-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 "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::Boundary& psiBf = psi.boundaryField();
183 
184  const fvMesh& mesh = psi.mesh();
185 
186  const dictionary& MULEScontrols = mesh.solverDict(psi.name());
187 
188  const label nLimiterIter
189  (
190  MULEScontrols.lookupOrDefault<label>("nLimiterIter", 3)
191  );
192 
193  const scalar smoothLimiter
194  (
195  MULEScontrols.lookupOrDefault<scalar>("smoothLimiter", 0)
196  );
197 
198  const scalar extremaCoeff
199  (
200  MULEScontrols.lookupOrDefault<scalar>("extremaCoeff", 0)
201  );
202 
203  const scalarField& psi0 = psi.oldTime();
204 
205  const labelUList& owner = mesh.owner();
206  const labelUList& neighb = mesh.neighbour();
207  tmp<volScalarField::Internal> tVsc = mesh.Vsc();
208  const scalarField& V = tVsc();
209 
210  const scalarField& phiBDIf = phiBD;
211  const surfaceScalarField::Boundary& phiBDBf =
212  phiBD.boundaryField();
213 
214  const scalarField& phiCorrIf = phiCorr;
215  const surfaceScalarField::Boundary& phiCorrBf =
216  phiCorr.boundaryField();
217 
219  (
220  IOobject
221  (
222  "lambda",
223  mesh.time().timeName(),
224  mesh,
225  IOobject::NO_READ,
226  IOobject::NO_WRITE,
227  false
228  ),
229  mesh,
230  dimless,
231  allLambda,
232  false // Use slices for the couples
233  );
234 
235  scalarField& lambdaIf = lambda;
236  surfaceScalarField::Boundary& lambdaBf = lambda.boundaryFieldRef();
237 
238  scalarField psiMaxn(psiIf.size(), psiMin);
239  scalarField psiMinn(psiIf.size(), psiMax);
240 
241  scalarField sumPhiBD(psiIf.size(), 0.0);
242 
243  scalarField sumPhip(psiIf.size(), 0.0);
244  scalarField mSumPhim(psiIf.size(), 0.0);
245 
246  forAll(phiCorrIf, facei)
247  {
248  const label own = owner[facei];
249  const label nei = neighb[facei];
250 
251  psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
252  psiMinn[own] = min(psiMinn[own], psiIf[nei]);
253 
254  psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
255  psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
256 
257  sumPhiBD[own] += phiBDIf[facei];
258  sumPhiBD[nei] -= phiBDIf[facei];
259 
260  const scalar phiCorrf = phiCorrIf[facei];
261 
262  if (phiCorrf > 0)
263  {
264  sumPhip[own] += phiCorrf;
265  mSumPhim[nei] += phiCorrf;
266  }
267  else
268  {
269  mSumPhim[own] -= phiCorrf;
270  sumPhip[nei] -= phiCorrf;
271  }
272  }
273 
274  forAll(phiCorrBf, patchi)
275  {
276  const fvPatchScalarField& psiPf = psiBf[patchi];
277  const scalarField& phiBDPf = phiBDBf[patchi];
278  const scalarField& phiCorrPf = phiCorrBf[patchi];
279 
280  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
281 
282  if (psiPf.coupled())
283  {
284  const scalarField psiPNf(psiPf.patchNeighbourField());
285 
286  forAll(phiCorrPf, pFacei)
287  {
288  const label pfCelli = pFaceCells[pFacei];
289 
290  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
291  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
292  }
293  }
294  else if (psiPf.fixesValue())
295  {
296  forAll(phiCorrPf, pFacei)
297  {
298  const label pfCelli = pFaceCells[pFacei];
299 
300  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
301  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
302  }
303  }
304  else
305  {
306  forAll(phiCorrPf, pFacei)
307  {
308  const label pfCelli = pFaceCells[pFacei];
309 
310  psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiMax);
311  psiMinn[pfCelli] = min(psiMinn[pfCelli], psiMin);
312  }
313  }
314 
315  forAll(phiCorrPf, pFacei)
316  {
317  const label pfCelli = pFaceCells[pFacei];
318 
319  sumPhiBD[pfCelli] += phiBDPf[pFacei];
320 
321  const scalar phiCorrf = phiCorrPf[pFacei];
322 
323  if (phiCorrf > 0)
324  {
325  sumPhip[pfCelli] += phiCorrf;
326  }
327  else
328  {
329  mSumPhim[pfCelli] -= phiCorrf;
330  }
331  }
332  }
333 
334  psiMaxn = min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
335  psiMinn = max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
336 
337  if (smoothLimiter > SMALL)
338  {
339  psiMaxn =
340  min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
341  psiMinn =
342  max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
343  }
344 
345  if (mesh.moving())
346  {
348 
349  psiMaxn =
350  V
351  *(
352  (rho.field()*rDeltaT - Sp.field())*psiMaxn
353  - Su.field()
354  )
355  - (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
356  + sumPhiBD;
357 
358  psiMinn =
359  V
360  *(
361  Su.field()
362  - (rho.field()*rDeltaT - Sp.field())*psiMinn
363  )
364  + (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
365  - sumPhiBD;
366  }
367  else
368  {
369  psiMaxn =
370  V
371  *(
372  (rho.field()*rDeltaT - Sp.field())*psiMaxn
373  - Su.field()
374  - (rho.oldTime().field()*rDeltaT)*psi0
375  )
376  + sumPhiBD;
377 
378  psiMinn =
379  V
380  *(
381  Su.field()
382  - (rho.field()*rDeltaT - Sp.field())*psiMinn
383  + (rho.oldTime().field()*rDeltaT)*psi0
384  )
385  - sumPhiBD;
386  }
387 
388  scalarField sumlPhip(psiIf.size());
389  scalarField mSumlPhim(psiIf.size());
390 
391  for (int j=0; j<nLimiterIter; j++)
392  {
393  sumlPhip = 0;
394  mSumlPhim = 0;
395 
396  forAll(lambdaIf, facei)
397  {
398  const label own = owner[facei];
399  const label nei = neighb[facei];
400 
401  scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
402 
403  if (lambdaPhiCorrf > 0)
404  {
405  sumlPhip[own] += lambdaPhiCorrf;
406  mSumlPhim[nei] += lambdaPhiCorrf;
407  }
408  else
409  {
410  mSumlPhim[own] -= lambdaPhiCorrf;
411  sumlPhip[nei] -= lambdaPhiCorrf;
412  }
413  }
414 
415  forAll(lambdaBf, patchi)
416  {
417  scalarField& lambdaPf = lambdaBf[patchi];
418  const scalarField& phiCorrfPf = phiCorrBf[patchi];
419 
420  const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
421 
422  forAll(lambdaPf, pFacei)
423  {
424  const label pfCelli = pFaceCells[pFacei];
425  const scalar lambdaPhiCorrf =
426  lambdaPf[pFacei]*phiCorrfPf[pFacei];
427 
428  if (lambdaPhiCorrf > 0)
429  {
430  sumlPhip[pfCelli] += lambdaPhiCorrf;
431  }
432  else
433  {
434  mSumlPhim[pfCelli] -= lambdaPhiCorrf;
435  }
436  }
437  }
438 
439  forAll(sumlPhip, celli)
440  {
441  sumlPhip[celli] =
442  max(min
443  (
444  (sumlPhip[celli] + psiMaxn[celli])
445  /(mSumPhim[celli] + ROOTVSMALL),
446  1.0), 0.0
447  );
448 
449  mSumlPhim[celli] =
450  max(min
451  (
452  (mSumlPhim[celli] + psiMinn[celli])
453  /(sumPhip[celli] + ROOTVSMALL),
454  1.0), 0.0
455  );
456  }
457 
458  const scalarField& lambdam = sumlPhip;
459  const scalarField& lambdap = mSumlPhim;
460 
461  forAll(lambdaIf, facei)
462  {
463  if (phiCorrIf[facei] > 0)
464  {
465  lambdaIf[facei] = min
466  (
467  lambdaIf[facei],
468  min(lambdap[owner[facei]], lambdam[neighb[facei]])
469  );
470  }
471  else
472  {
473  lambdaIf[facei] = min
474  (
475  lambdaIf[facei],
476  min(lambdam[owner[facei]], lambdap[neighb[facei]])
477  );
478  }
479  }
480 
481  forAll(lambdaBf, patchi)
482  {
483  fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
484  const scalarField& phiCorrfPf = phiCorrBf[patchi];
485  const fvPatchScalarField& psiPf = psiBf[patchi];
486 
487  if (isA<wedgeFvPatch>(mesh.boundary()[patchi]))
488  {
489  lambdaPf = 0;
490  }
491  else if (psiPf.coupled())
492  {
493  const labelList& pFaceCells =
494  mesh.boundary()[patchi].faceCells();
495 
496  forAll(lambdaPf, pFacei)
497  {
498  const label pfCelli = pFaceCells[pFacei];
499 
500  if (phiCorrfPf[pFacei] > 0)
501  {
502  lambdaPf[pFacei] =
503  min(lambdaPf[pFacei], lambdap[pfCelli]);
504  }
505  else
506  {
507  lambdaPf[pFacei] =
508  min(lambdaPf[pFacei], lambdam[pfCelli]);
509  }
510  }
511  }
512  }
513 
514  syncTools::syncFaceList(mesh, allLambda, minEqOp<scalar>());
515  }
516 }
517 
518 
519 template<class RdeltaTType, class RhoType, class SpType, class SuType>
521 (
522  const RdeltaTType& rDeltaT,
523  const RhoType& rho,
524  const volScalarField& psi,
525  const surfaceScalarField& phi,
526  surfaceScalarField& phiPsi,
527  const SpType& Sp,
528  const SuType& Su,
529  const scalar psiMax,
530  const scalar psiMin,
531  const bool returnCorr
532 )
533 {
534  const fvMesh& mesh = psi.mesh();
535 
536  surfaceScalarField phiBD(upwind<scalar>(psi.mesh(), phi).flux(psi));
537 
538  surfaceScalarField::Boundary& phiBDBf = phiBD.boundaryFieldRef();
539  const surfaceScalarField::Boundary& phiPsiBf = phiPsi.boundaryField();
540 
541  forAll(phiBDBf, patchi)
542  {
543  fvsPatchScalarField& phiBDPf = phiBDBf[patchi];
544 
545  if (!phiBDPf.coupled())
546  {
547  phiBDPf = phiPsiBf[patchi];
548  }
549  }
550 
551  surfaceScalarField& phiCorr = phiPsi;
552  phiCorr -= phiBD;
553 
554  scalarField allLambda(mesh.nFaces(), 1.0);
555 
557  (
558  IOobject
559  (
560  "lambda",
561  mesh.time().timeName(),
562  mesh,
563  IOobject::NO_READ,
564  IOobject::NO_WRITE,
565  false
566  ),
567  mesh,
568  dimless,
569  allLambda,
570  false // Use slices for the couples
571  );
572 
573  limiter
574  (
575  allLambda,
576  rDeltaT,
577  rho,
578  psi,
579  phiBD,
580  phiCorr,
581  Sp,
582  Su,
583  psiMax,
584  psiMin
585  );
586 
587  if (returnCorr)
588  {
589  phiCorr *= lambda;
590  }
591  else
592  {
593  phiPsi = phiBD + lambda*phiCorr;
594  }
595 }
596 
597 
598 template<class SurfaceScalarFieldList>
599 void Foam::MULES::limitSum(SurfaceScalarFieldList& phiPsiCorrs)
600 {
601  {
602  UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
603  forAll(phiPsiCorrs, phasei)
604  {
605  phiPsiCorrsInternal.set(phasei, &phiPsiCorrs[phasei]);
606  }
607 
608  limitSum(phiPsiCorrsInternal);
609  }
610 
611  const surfaceScalarField::Boundary& bfld =
612  phiPsiCorrs[0].boundaryField();
613 
614  forAll(bfld, patchi)
615  {
616  if (bfld[patchi].coupled())
617  {
618  UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
619  forAll(phiPsiCorrs, phasei)
620  {
621  phiPsiCorrsPatch.set
622  (
623  phasei,
624  &phiPsiCorrs[phasei].boundaryFieldRef()[patchi]
625  );
626  }
627 
628  limitSum(phiPsiCorrsPatch);
629  }
630  }
631 }
632 
633 
634 // ************************************************************************* //
surfaceScalarField & phi
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
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)
#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:291
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
bool moving() const
Is mesh moving.
Definition: polyMesh.H:496
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:253
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: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
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:51
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)
virtual bool coupled() const
Return true if this patch field is coupled.
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