PhiScheme.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 "volFields.H"
27 #include "surfaceFields.H"
28 #include "fvcGrad.H"
29 #include "coupledFvPatchFields.H"
30 #include "surfaceInterpolate.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 template<class Type, class PhiLimiter>
37 (
39 ) const
40 {
41  const fvMesh& mesh = this->mesh();
42 
44  (
46  (
47  IOobject
48  (
49  "PhiLimiter",
50  mesh.time().timeName(),
51  mesh
52  ),
53  mesh,
54  dimless
55  )
56  );
57  surfaceScalarField& Limiter = tLimiter.ref();
58 
59  const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights();
60 
61  const surfaceVectorField& Sf = mesh.Sf();
62  const surfaceScalarField& magSf = mesh.magSf();
63 
64  const labelUList& owner = mesh.owner();
65  const labelUList& neighbour = mesh.neighbour();
66 
67  tmp<surfaceScalarField> tUflux = this->faceFlux_;
68 
69  if (this->faceFlux_.dimensions() == dimDensity*dimVelocity*dimArea)
70  {
71  const volScalarField& rho =
72  phi.db().objectRegistry::template lookupObject<volScalarField>
73  ("rho");
74 
75  tUflux = this->faceFlux_/fvc::interpolate(rho);
76  }
77  else if (this->faceFlux_.dimensions() != dimVelocity*dimArea)
78  {
80  << "dimensions of faceFlux are not correct"
81  << exit(FatalError);
82  }
83 
84  const surfaceScalarField& Uflux = tUflux();
85 
86  scalarField& pLimiter = Limiter.primitiveFieldRef();
87 
88  forAll(pLimiter, face)
89  {
90  pLimiter[face] = PhiLimiter::limiter
91  (
92  CDweights[face],
93  Uflux[face],
94  phi[owner[face]],
95  phi[neighbour[face]],
96  Sf[face],
97  magSf[face]
98  );
99  }
100 
101 
102  surfaceScalarField::Boundary& bLimiter =
103  Limiter.boundaryFieldRef();
104 
105  forAll(bLimiter, patchi)
106  {
107  scalarField& pLimiter = bLimiter[patchi];
108 
109  if (bLimiter[patchi].coupled())
110  {
111  const scalarField& pCDweights = CDweights.boundaryField()[patchi];
112  const vectorField& pSf = Sf.boundaryField()[patchi];
113  const scalarField& pmagSf = magSf.boundaryField()[patchi];
114  const scalarField& pFaceFlux = Uflux.boundaryField()[patchi];
115 
116  const Field<Type> pphiP
117  (
118  phi.boundaryField()[patchi].patchInternalField()
119  );
120  const Field<Type> pphiN
121  (
122  phi.boundaryField()[patchi].patchNeighbourField()
123  );
124 
125  forAll(pLimiter, face)
126  {
127  pLimiter[face] = PhiLimiter::limiter
128  (
129  pCDweights[face],
130  pFaceFlux[face],
131  pphiP[face],
132  pphiN[face],
133  pSf[face],
134  pmagSf[face]
135  );
136  }
137  }
138  else
139  {
140  pLimiter = 1.0;
141  }
142  }
143 
144  return tLimiter;
145 }
146 
147 
148 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const surfaceVectorField & Sf() const
Return cell face area vectors.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Generic GeometricField class.
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
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Definition: PhiScheme.C:37
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:288
dynamicFvMesh & mesh
Calculate the gradient of the given field.
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)
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
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
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
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:361
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
const dimensionSet dimVelocity