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-2021 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  "PhiLimiter",
48  mesh,
49  dimless
50  )
51  );
52  surfaceScalarField& Limiter = tLimiter.ref();
53 
54  const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights();
55 
56  const surfaceVectorField& Sf = mesh.Sf();
57  const surfaceScalarField& magSf = mesh.magSf();
58 
59  const labelUList& owner = mesh.owner();
60  const labelUList& neighbour = mesh.neighbour();
61 
62  tmp<surfaceScalarField> tUflux = this->faceFlux_;
63 
64  if (this->faceFlux_.dimensions() == dimMassFlux)
65  {
66  const volScalarField& rho =
67  phi.db().objectRegistry::template lookupObject<volScalarField>
68  ("rho");
69 
70  tUflux = this->faceFlux_/fvc::interpolate(rho);
71  }
72  else if (this->faceFlux_.dimensions() != dimFlux)
73  {
75  << "dimensions of faceFlux are not correct"
76  << exit(FatalError);
77  }
78 
79  const surfaceScalarField& Uflux = tUflux();
80 
81  scalarField& pLimiter = Limiter.primitiveFieldRef();
82 
83  forAll(pLimiter, face)
84  {
85  pLimiter[face] = PhiLimiter::limiter
86  (
87  CDweights[face],
88  Uflux[face],
89  phi[owner[face]],
90  phi[neighbour[face]],
91  Sf[face],
92  magSf[face]
93  );
94  }
95 
96 
97  surfaceScalarField::Boundary& bLimiter =
98  Limiter.boundaryFieldRef();
99 
100  forAll(bLimiter, patchi)
101  {
102  scalarField& pLimiter = bLimiter[patchi];
103 
104  if (bLimiter[patchi].coupled())
105  {
106  const scalarField& pCDweights = CDweights.boundaryField()[patchi];
107  const vectorField& pSf = Sf.boundaryField()[patchi];
108  const scalarField& pmagSf = magSf.boundaryField()[patchi];
109  const scalarField& pFaceFlux = Uflux.boundaryField()[patchi];
110 
111  const Field<Type> pphiP
112  (
113  phi.boundaryField()[patchi].patchInternalField()
114  );
115  const Field<Type> pphiN
116  (
117  phi.boundaryField()[patchi].patchNeighbourField()
118  );
119 
120  forAll(pLimiter, face)
121  {
122  pLimiter[face] = PhiLimiter::limiter
123  (
124  pCDweights[face],
125  pFaceFlux[face],
126  pphiP[face],
127  pphiN[face],
128  pSf[face],
129  pmagSf[face]
130  );
131  }
132  }
133  else
134  {
135  pLimiter = 1.0;
136  }
137  }
138 
139  return tLimiter;
140 }
141 
142 
143 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:323
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Generic GeometricField class.
const dimensionSet dimless
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:284
dynamicFvMesh & mesh
Calculate the gradient of the given field.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
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)
const dimensionSet dimFlux
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:60
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
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 labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:278
const dimensionSet dimMassFlux
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
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:312