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-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 "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  (
44  "tmp<surfaceScalarField> interpolate(const RhoType& rho)"
45  );
46  return tmp<surfaceScalarField>(NULL);
47  }
48 
49  template<>
51  {
52  return fvc::interpolate(rho);
53  }
54 }
55 }
56 
57 
58 template<class RhoType, class SpType, class SuType>
60 (
61  const RhoType& rho,
62  volScalarField& psi,
63  const surfaceScalarField& phi,
64  surfaceScalarField& phiPsi,
65  const SpType& Sp,
66  const SuType& Su,
67  const scalar psiMax,
68  const scalar psiMin
69 )
70 {
71  const fvMesh& mesh = psi.mesh();
72 
73  const dictionary& MULEScontrols = mesh.solverDict(psi.name());
74 
75  label maxIter
76  (
77  readLabel(MULEScontrols.lookup("maxIter"))
78  );
79 
80  scalar maxUnboundedness
81  (
82  readScalar(MULEScontrols.lookup("maxUnboundedness"))
83  );
84 
85  scalar CoCoeff
86  (
87  readScalar(MULEScontrols.lookup("CoCoeff"))
88  );
89 
90  scalarField allCoLambda(mesh.nFaces());
91 
92  {
94  (
95  IOobject
96  (
97  "CoLambda",
98  mesh.time().timeName(),
99  mesh,
102  false
103  ),
104  mesh,
105  dimless,
106  allCoLambda,
107  false // Use slices for the couples
108  );
109 
111  {
113  mesh.time().deltaT()*mesh.surfaceInterpolation::deltaCoeffs()
114  *mag(phi/interpolate(rho))/mesh.magSf();
115 
116  CoLambda == 1.0/max(CoCoeff*Cof, scalar(1));
117  }
118  else
119  {
121  mesh.time().deltaT()*mesh.surfaceInterpolation::deltaCoeffs()
122  *mag(phi)/mesh.magSf();
123 
124  CoLambda == 1.0/max(CoCoeff*Cof, scalar(1));
125  }
126  }
127 
128  scalarField allLambda(allCoLambda);
129  //scalarField allLambda(mesh.nFaces(), 1.0);
130 
132  (
133  IOobject
134  (
135  "lambda",
136  mesh.time().timeName(),
137  mesh,
140  false
141  ),
142  mesh,
143  dimless,
144  allLambda,
145  false // Use slices for the couples
146  );
147 
148  linear<scalar> CDs(mesh);
149  upwind<scalar> UDs(mesh, phi);
150  //fv::uncorrectedSnGrad<scalar> snGrads(mesh);
151 
152  fvScalarMatrix psiConvectionDiffusion
153  (
154  fvm::ddt(rho, psi)
155  + fv::gaussConvectionScheme<scalar>(mesh, phi, UDs).fvmDiv(phi, psi)
156  //- fv::gaussLaplacianScheme<scalar, scalar>(mesh, CDs, snGrads)
157  //.fvmLaplacian(Dpsif, psi)
158  - fvm::Sp(Sp, psi)
159  - Su
160  );
161 
162  surfaceScalarField phiBD(psiConvectionDiffusion.flux());
163 
164  surfaceScalarField& phiCorr = phiPsi;
165  phiCorr -= phiBD;
166 
167  for (label i=0; i<maxIter; i++)
168  {
169  if (i != 0 && i < 4)
170  {
171  allLambda = allCoLambda;
172  }
173 
174  limiter
175  (
176  allLambda,
177  1.0/mesh.time().deltaTValue(),
178  rho,
179  psi,
180  phiBD,
181  phiCorr,
182  Sp,
183  Su,
184  psiMax,
185  psiMin
186  );
187 
188  solve
189  (
190  psiConvectionDiffusion + fvc::div(lambda*phiCorr),
191  MULEScontrols
192  );
193 
194  scalar maxPsiM1 = gMax(psi.internalField()) - 1.0;
195  scalar minPsi = gMin(psi.internalField());
196 
197  scalar unboundedness = max(max(maxPsiM1, 0.0), -min(minPsi, 0.0));
198 
199  if (unboundedness < maxUnboundedness)
200  {
201  break;
202  }
203  else
204  {
205  Info<< "MULES: max(" << psi.name() << " - 1) = " << maxPsiM1
206  << " min(" << psi.name() << ") = " << minPsi << endl;
207 
208  phiBD = psiConvectionDiffusion.flux();
209 
210  /*
211  word gammaScheme("div(phi,gamma)");
212  word gammarScheme("div(phirb,gamma)");
213 
214  const surfaceScalarField& phir =
215  mesh.lookupObject<surfaceScalarField>("phir");
216 
217  phiCorr =
218  fvc::flux
219  (
220  phi,
221  psi,
222  gammaScheme
223  )
224  + fvc::flux
225  (
226  -fvc::flux(-phir, scalar(1) - psi, gammarScheme),
227  psi,
228  gammarScheme
229  )
230  - phiBD;
231  */
232  }
233  }
234 
235  phiPsi = psiConvectionDiffusion.flux() + lambda*phiCorr;
236 }
237 
238 
239 // ************************************************************************* //
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label nFaces() const
Calculate the matrix for implicit and explicit sources.
dimensioned< scalar > mag(const dimensioned< Type > &)
Central-differencing interpolation scheme class.
Definition: linear.H:50
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:883
IMULES: Multidimensional universal limiter for implicit solution.
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
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
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:741
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:365
Calulate the matrix for the first temporal derivative.
Surface Interpolation.
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:68
messageStream Info
tmp< surfaceScalarField > interpolate(const RhoType &rho)
const Mesh & mesh() const
Return mesh.
dynamicFvMesh & mesh
const dimensionSet dimDensity
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Namespace for OpenFOAM.
label readLabel(Istream &is)
Definition: label.H:64
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
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
Basic second-order convection using face-gradients and Gauss&#39; theorem.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
const volScalarField & psi
Definition: createFields.H:24
const dimensionSet & dimensions() const
Return dimensions.
const word & name() const
Return name.
Definition: IOobject.H:260
Type gMin(const FieldField< Field, Type > &f)
solverPerformance solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Calculate the divergence of the given field.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
Upwind differencing scheme class.
Definition: upwind.H:51
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
dimensionedScalar deltaT() const
Return time step.
Definition: TimeState.C:79
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:356
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)
const dimensionSet dimVelocity
void implicitSolve(const RhoType &rho, volScalarField &gamma, const surfaceScalarField &phi, surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
Type gMax(const FieldField< Field, Type > &f)
A class for managing temporary objects.
Definition: PtrList.H:118