IMULESTemplates.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-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 "IMULES.H"
27 #include "gaussConvectionScheme.H"
28 #include "surfaceInterpolate.H"
29 #include "fvmDdt.H"
30 #include "fvmSup.H"
31 #include "fvcDiv.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace MULES
38 {
39  template<class RhoType>
40  inline tmp<surfaceScalarField> interpolate(const RhoType& rho)
41  {
43  return tmp<surfaceScalarField>(NULL);
44  }
45 
46  template<>
48  {
49  return fvc::interpolate(rho);
50  }
51 }
52 }
53 
54 
55 template<class RhoType, class SpType, class SuType>
57 (
58  const RhoType& rho,
59  volScalarField& psi,
60  const surfaceScalarField& phi,
61  surfaceScalarField& phiPsi,
62  const SpType& Sp,
63  const SuType& Su,
64  const scalar psiMax,
65  const scalar psiMin
66 )
67 {
68  const fvMesh& mesh = psi.mesh();
69 
70  const dictionary& MULEScontrols = mesh.solverDict(psi.name());
71 
72  label maxIter
73  (
74  readLabel(MULEScontrols.lookup("maxIter"))
75  );
76 
77  scalar maxUnboundedness
78  (
79  readScalar(MULEScontrols.lookup("maxUnboundedness"))
80  );
81 
82  scalar CoCoeff
83  (
84  readScalar(MULEScontrols.lookup("CoCoeff"))
85  );
86 
87  scalarField allCoLambda(mesh.nFaces());
88 
89  {
91  (
92  IOobject
93  (
94  "CoLambda",
95  mesh.time().timeName(),
96  mesh,
99  false
100  ),
101  mesh,
102  dimless,
103  allCoLambda,
104  false // Use slices for the couples
105  );
106 
108  {
110  mesh.time().deltaT()*mesh.surfaceInterpolation::deltaCoeffs()
111  *mag(phi/interpolate(rho))/mesh.magSf();
112 
113  CoLambda == 1.0/max(CoCoeff*Cof, scalar(1));
114  }
115  else
116  {
118  mesh.time().deltaT()*mesh.surfaceInterpolation::deltaCoeffs()
119  *mag(phi)/mesh.magSf();
120 
121  CoLambda == 1.0/max(CoCoeff*Cof, scalar(1));
122  }
123  }
124 
125  scalarField allLambda(allCoLambda);
126  //scalarField allLambda(mesh.nFaces(), 1.0);
127 
129  (
130  IOobject
131  (
132  "lambda",
133  mesh.time().timeName(),
134  mesh,
137  false
138  ),
139  mesh,
140  dimless,
141  allLambda,
142  false // Use slices for the couples
143  );
144 
145  linear<scalar> CDs(mesh);
146  upwind<scalar> UDs(mesh, phi);
147  //fv::uncorrectedSnGrad<scalar> snGrads(mesh);
148 
149  fvScalarMatrix psiConvectionDiffusion
150  (
151  fvm::ddt(rho, psi)
152  + fv::gaussConvectionScheme<scalar>(mesh, phi, UDs).fvmDiv(phi, psi)
153  //- fv::gaussLaplacianScheme<scalar, scalar>(mesh, CDs, snGrads)
154  //.fvmLaplacian(Dpsif, psi)
155  - fvm::Sp(Sp, psi)
156  - Su
157  );
158 
159  surfaceScalarField phiBD(psiConvectionDiffusion.flux());
160 
161  surfaceScalarField& phiCorr = phiPsi;
162  phiCorr -= phiBD;
163 
164  for (label i=0; i<maxIter; i++)
165  {
166  if (i != 0 && i < 4)
167  {
168  allLambda = allCoLambda;
169  }
170 
171  limiter
172  (
173  allLambda,
174  1.0/mesh.time().deltaTValue(),
175  rho,
176  psi,
177  phiBD,
178  phiCorr,
179  Sp,
180  Su,
181  psiMax,
182  psiMin
183  );
184 
185  solve
186  (
187  psiConvectionDiffusion + fvc::div(lambda*phiCorr),
188  MULEScontrols
189  );
190 
191  scalar maxPsiM1 = gMax(psi.primitiveField()) - 1.0;
192  scalar minPsi = gMin(psi.primitiveField());
193 
194  scalar unboundedness = max(max(maxPsiM1, 0.0), -min(minPsi, 0.0));
195 
196  if (unboundedness < maxUnboundedness)
197  {
198  break;
199  }
200  else
201  {
202  Info<< "MULES: max(" << psi.name() << " - 1) = " << maxPsiM1
203  << " min(" << psi.name() << ") = " << minPsi << endl;
204 
205  phiBD = psiConvectionDiffusion.flux();
206 
207  /*
208  word gammaScheme("div(phi,gamma)");
209  word gammarScheme("div(phirb,gamma)");
210 
211  const surfaceScalarField& phir =
212  mesh.lookupObject<surfaceScalarField>("phir");
213 
214  phiCorr =
215  fvc::flux
216  (
217  phi,
218  psi,
219  gammaScheme
220  )
221  + fvc::flux
222  (
223  -fvc::flux(-phir, scalar(1) - psi, gammarScheme),
224  psi,
225  gammarScheme
226  )
227  - phiBD;
228  */
229  }
230  }
231 
232  phiPsi = psiConvectionDiffusion.flux() + lambda*phiCorr;
233 }
234 
235 
236 // ************************************************************************* //
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
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)
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 > &)
Type gMin(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Central-differencing interpolation scheme class.
Definition: linear.H:50
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:864
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
tmp< surfaceScalarField > interpolate(const RhoType &rho)
dynamicFvMesh & mesh
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calulate the matrix for the first temporal derivative.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:71
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
label readLabel(Istream &is)
Definition: label.H:64
void implicitSolve(const RhoType &rho, volScalarField &gamma, const surfaceScalarField &phi, surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
const dimensionSet & dimensions() const
Return dimensions.
Calculate the divergence of the given field.
Type gMax(const FieldField< Field, Type > &f)
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
const dimensionSet dimDensity
label nFaces() const
Basic second-order convection using face-gradients and Gauss&#39; theorem.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:53
const Mesh & mesh() const
Return mesh.
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"))
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const volScalarField & psi
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:342
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 ...
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const word & name() const
Return name.
Definition: IOobject.H:260
Calculate the matrix for implicit and explicit sources.
IMULES: Multidimensional universal limiter for implicit solution.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
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
const dimensionSet dimVelocity