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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "fvcGrad.H"
32 #include "coupledFvPatchFields.H"
33 
34 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
35 
36 template<class Type, class Limiter, template<class> class LimitFunc>
38 (
39  const GeometricField<Type, fvPatchField, volMesh>& phi,
40  surfaceScalarField& limiterField
41 ) const
42 {
43  typedef GeometricField<typename Limiter::phiType, fvPatchField, volMesh>
44  VolFieldType;
45 
46  typedef GeometricField<typename Limiter::gradPhiType, fvPatchField, volMesh>
47  GradVolFieldType;
48 
49  const fvMesh& mesh = this->mesh();
50 
51  tmp<VolFieldType> tlPhi = LimitFunc<Type>()(phi);
52  const VolFieldType& lPhi = tlPhi();
53 
54  tmp<GradVolFieldType> tgradc(fvc::grad(lPhi));
55  const GradVolFieldType& gradc = tgradc();
56 
57  const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights();
58 
59  const labelUList& owner = mesh.owner();
60  const labelUList& neighbour = mesh.neighbour();
61 
62  const vectorField& C = mesh.C();
63 
64  scalarField& pLim = limiterField.primitiveFieldRef();
65 
66  forAll(pLim, face)
67  {
68  label own = owner[face];
69  label nei = neighbour[face];
70 
71  pLim[face] = Limiter::limiter
72  (
73  CDweights[face],
74  this->faceFlux_[face],
75  lPhi[own],
76  lPhi[nei],
77  gradc[own],
78  gradc[nei],
79  C[nei] - C[own]
80  );
81  }
82 
83  surfaceScalarField::Boundary& bLim = limiterField.boundaryFieldRef();
84 
85  forAll(bLim, patchi)
86  {
87  scalarField& pLim = bLim[patchi];
88 
89  if (bLim[patchi].coupled())
90  {
91  const scalarField& pCDweights = CDweights.boundaryField()[patchi];
92  const scalarField& pFaceFlux =
93  this->faceFlux_.boundaryField()[patchi];
94 
95  const Field<typename Limiter::phiType> plPhiP
96  (
97  lPhi.boundaryField()[patchi].patchInternalField()
98  );
99  const Field<typename Limiter::phiType> plPhiN
100  (
101  lPhi.boundaryField()[patchi].patchNeighbourField()
102  );
103  const Field<typename Limiter::gradPhiType> pGradcP
104  (
105  gradc.boundaryField()[patchi].patchInternalField()
106  );
107  const Field<typename Limiter::gradPhiType> pGradcN
108  (
109  gradc.boundaryField()[patchi].patchNeighbourField()
110  );
111 
112  // Build the d-vectors
113  vectorField pd(CDweights.boundaryField()[patchi].patch().delta());
114 
115  forAll(pLim, face)
116  {
117  pLim[face] = Limiter::limiter
118  (
119  pCDweights[face],
120  pFaceFlux[face],
121  plPhiP[face],
122  plPhiN[face],
123  pGradcP[face],
124  pGradcN[face],
125  pd[face]
126  );
127  }
128  }
129  else
130  {
131  pLim = 1.0;
132  }
133  }
134 
135  limiterField.setOriented();
136 }
138 
139 // * * * * * * * * * * * * Public Member Functions * * * * * * * * * * * * //
140 
141 template<class Type, class Limiter, template<class> class LimitFunc>
144 (
146 ) const
147 {
148  const fvMesh& mesh = this->mesh();
149 
150  const word limiterFieldName(type() + "Limiter(" + phi.name() + ')');
151 
152  if (this->mesh().cache("limiter"))
153  {
154  auto* fldptr = mesh.getObjectPtr<surfaceScalarField>(limiterFieldName);
155 
156  if (!fldptr)
157  {
158  fldptr = new surfaceScalarField
159  (
160  IOobject
161  (
162  limiterFieldName,
163  mesh.time().timeName(),
164  mesh,
165  IOobject::NO_READ,
166  IOobject::NO_WRITE,
167  IOobject::REGISTER
168  ),
169  mesh,
170  dimless
171  );
172 
173  mesh.objectRegistry::store(fldptr);
174  }
175  auto& limiterField = *fldptr;
176 
177  calcLimiter(phi, limiterField);
178 
180  (
181  limiterFieldName,
182  limiterField
183  );
184  }
185  else
186  {
187  tmp<surfaceScalarField> tlimiterField
188  (
190  (
191  IOobject
192  (
193  limiterFieldName,
194  mesh.time().timeName(),
195  mesh
196  ),
197  mesh,
198  dimless
199  )
200  );
201 
202  calcLimiter(phi, tlimiterField.ref());
203 
204  return tlimiterField;
205  }
206 }
207 
208 
209 // ************************************************************************* //
Foam::surfaceFields.
type
Types of root.
Definition: Roots.H:52
Class to create NVD/TVD limited weighting-factors.
Definition: LimitedScheme.H:64
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
const dimensionSet dimless
Dimensionless.
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:31
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
Calculate the gradient of the given field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
volScalarField & C
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field< vector > vectorField
Specialisation of Field<T> for vector.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
bool coupled
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51