PhiScheme.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "fvcGrad.H"
31 #include "coupledFvPatchFields.H"
32 #include "surfaceInterpolate.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 template<class Type, class PhiLimiter>
39 (
41 ) const
42 {
43  const fvMesh& mesh = this->mesh();
44 
46  (
48  (
49  IOobject
50  (
51  "PhiLimiter",
52  mesh.time().timeName(),
53  mesh
54  ),
55  mesh,
56  dimless
57  )
58  );
59  surfaceScalarField& Limiter = tLimiter.ref();
60 
61  const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights();
62 
63  const surfaceVectorField& Sf = mesh.Sf();
64  const surfaceScalarField& magSf = mesh.magSf();
65 
66  const labelUList& owner = mesh.owner();
67  const labelUList& neighbour = mesh.neighbour();
68 
69  tmp<surfaceScalarField> tUflux = this->faceFlux_;
70 
71  if (this->faceFlux_.dimensions() == dimDensity*dimVelocity*dimArea)
72  {
73  const volScalarField& rho =
74  phi.db().objectRegistry::template lookupObject<volScalarField>
75  ("rho");
76 
77  tUflux = this->faceFlux_/fvc::interpolate(rho);
78  }
79  else if (this->faceFlux_.dimensions() != dimVelocity*dimArea)
80  {
82  << "dimensions of faceFlux are not correct"
83  << exit(FatalError);
84  }
85 
86  const surfaceScalarField& Uflux = tUflux();
87 
88  scalarField& pLimiter = Limiter.primitiveFieldRef();
89 
90  forAll(pLimiter, face)
91  {
92  pLimiter[face] = PhiLimiter::limiter
93  (
94  CDweights[face],
95  Uflux[face],
96  phi[owner[face]],
97  phi[neighbour[face]],
98  Sf[face],
99  magSf[face]
100  );
101  }
102 
103 
104  surfaceScalarField::Boundary& bLimiter =
105  Limiter.boundaryFieldRef();
106 
107  forAll(bLimiter, patchi)
108  {
109  scalarField& pLimiter = bLimiter[patchi];
110 
111  if (bLimiter[patchi].coupled())
112  {
113  const scalarField& pCDweights = CDweights.boundaryField()[patchi];
114  const vectorField& pSf = Sf.boundaryField()[patchi];
115  const scalarField& pmagSf = magSf.boundaryField()[patchi];
116  const scalarField& pFaceFlux = Uflux.boundaryField()[patchi];
117 
118  const Field<Type> pphiP
119  (
120  phi.boundaryField()[patchi].patchInternalField()
121  );
122  const Field<Type> pphiN
123  (
124  phi.boundaryField()[patchi].patchNeighbourField()
125  );
126 
127  forAll(pLimiter, face)
128  {
129  pLimiter[face] = PhiLimiter::limiter
130  (
131  pCDweights[face],
132  pFaceFlux[face],
133  pphiP[face],
134  pphiN[face],
135  pSf[face],
136  pmagSf[face]
137  );
138  }
139  }
140  else
141  {
142  pLimiter = 1.0;
143  }
144  }
145 
146  return tLimiter;
147 }
148 
149 
150 // ************************************************************************* //
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
surfaceFields.H
Foam::surfaceFields.
rho
rho
Definition: readInitialConditions.H:88
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
coupledFvPatchFields.H
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::Field< scalar >
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::PhiScheme::limiter
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Definition: PhiScheme.C:39
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
fvcGrad.H
Calculate the gradient of the given field.
Foam::UList< label >
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::limiter
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:37
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189