LimitedScheme.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 
33 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
34 
35 template<class Type, class Limiter, template<class> class LimitFunc>
37 (
38  const GeometricField<Type, fvPatchField, volMesh>& phi,
39  surfaceScalarField& limiterField
40 ) const
41 {
42  typedef GeometricField<typename Limiter::phiType, fvPatchField, volMesh>
43  VolFieldType;
44 
45  typedef GeometricField<typename Limiter::gradPhiType, fvPatchField, volMesh>
46  GradVolFieldType;
47 
48  const fvMesh& mesh = this->mesh();
49 
50  tmp<VolFieldType> tlPhi = LimitFunc<Type>()(phi);
51  const VolFieldType& lPhi = tlPhi();
52 
53  tmp<GradVolFieldType> tgradc(fvc::grad(lPhi));
54  const GradVolFieldType& gradc = tgradc();
55 
56  const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights();
57 
58  const labelUList& owner = mesh.owner();
59  const labelUList& neighbour = mesh.neighbour();
60 
61  const vectorField& C = mesh.C();
62 
63  scalarField& pLim = limiterField.primitiveFieldRef();
64 
65  forAll(pLim, face)
66  {
67  label own = owner[face];
68  label nei = neighbour[face];
69 
70  pLim[face] = Limiter::limiter
71  (
72  CDweights[face],
73  this->faceFlux_[face],
74  lPhi[own],
75  lPhi[nei],
76  gradc[own],
77  gradc[nei],
78  C[nei] - C[own]
79  );
80  }
81 
82  surfaceScalarField::Boundary& bLim = limiterField.boundaryFieldRef();
83 
84  forAll(bLim, patchi)
85  {
86  scalarField& pLim = bLim[patchi];
87 
88  if (bLim[patchi].coupled())
89  {
90  const scalarField& pCDweights = CDweights.boundaryField()[patchi];
91  const scalarField& pFaceFlux =
92  this->faceFlux_.boundaryField()[patchi];
93 
94  const Field<typename Limiter::phiType> plPhiP
95  (
96  lPhi.boundaryField()[patchi].patchInternalField()
97  );
98  const Field<typename Limiter::phiType> plPhiN
99  (
100  lPhi.boundaryField()[patchi].patchNeighbourField()
101  );
102  const Field<typename Limiter::gradPhiType> pGradcP
103  (
104  gradc.boundaryField()[patchi].patchInternalField()
105  );
106  const Field<typename Limiter::gradPhiType> pGradcN
107  (
108  gradc.boundaryField()[patchi].patchNeighbourField()
109  );
110 
111  // Build the d-vectors
112  vectorField pd(CDweights.boundaryField()[patchi].patch().delta());
113 
114  forAll(pLim, face)
115  {
116  pLim[face] = Limiter::limiter
117  (
118  pCDweights[face],
119  pFaceFlux[face],
120  plPhiP[face],
121  plPhiN[face],
122  pGradcP[face],
123  pGradcN[face],
124  pd[face]
125  );
126  }
127  }
128  else
129  {
130  pLim = 1.0;
131  }
132  }
133 
134  limiterField.setOriented();
135 }
136 
137 
138 // * * * * * * * * * * * * Public Member Functions * * * * * * * * * * * * //
139 
140 template<class Type, class Limiter, template<class> class LimitFunc>
143 (
145 ) const
146 {
147  const fvMesh& mesh = this->mesh();
148 
149  const word limiterFieldName(type() + "Limiter(" + phi.name() + ')');
150 
151  if (this->mesh().cache("limiter"))
152  {
153  if (!mesh.foundObject<surfaceScalarField>(limiterFieldName))
154  {
155  surfaceScalarField* limiterField
156  (
158  (
159  IOobject
160  (
161  limiterFieldName,
162  mesh.time().timeName(),
163  mesh,
164  IOobject::NO_READ,
165  IOobject::NO_WRITE
166  ),
167  mesh,
168  dimless
169  )
170  );
171 
172  mesh.objectRegistry::store(limiterField);
173  }
174 
175  surfaceScalarField& limiterField =
176  mesh.lookupObjectRef<surfaceScalarField>(limiterFieldName);
177 
178  calcLimiter(phi, limiterField);
179 
180  return limiterField;
181  }
182  else
183  {
184  tmp<surfaceScalarField> tlimiterField
185  (
187  (
188  IOobject
189  (
190  limiterFieldName,
191  mesh.time().timeName(),
192  mesh
193  ),
194  mesh,
195  dimless
196  )
197  );
198 
199  calcLimiter(phi, tlimiterField.ref());
200 
201  return tlimiterField;
202  }
203 }
204 
205 
206 // ************************************************************************* //
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::LimitedScheme::limiter
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Definition: LimitedScheme.C:143
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
surfaceFields.H
Foam::surfaceFields.
Foam::LimitedScheme
Class to create NVD/TVD limited weighting-factors.
Definition: LimitedScheme.H:68
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
C
volScalarField & C
Definition: readThermalProperties.H:102
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
coupledFvPatchFields.H
Foam::label
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:62
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
fvcGrad.H
Calculate the gradient of the given field.
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:79
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::limiter
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:37