LimitedScheme.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-2022 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 
31 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
32 
33 template<class Type, class Limiter, template<class> class LimitFunc>
35 (
36  const GeometricField<Type, fvPatchField, volMesh>& phi,
37  surfaceScalarField& limiterField
38 ) const
39 {
40  const fvMesh& mesh = this->mesh();
41 
42  tmp<GeometricField<typename Limiter::phiType, fvPatchField, volMesh>>
43  tlPhi = LimitFunc<Type>()(phi);
44 
45  const GeometricField<typename Limiter::phiType, fvPatchField, volMesh>&
46  lPhi = tlPhi();
47 
48  tmp<GeometricField<typename Limiter::gradPhiType, fvPatchField, volMesh>>
49  tgradc(fvc::grad(lPhi));
50  const GeometricField<typename Limiter::gradPhiType, fvPatchField, volMesh>&
51  gradc = tgradc();
52 
53  const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights();
54 
55  const labelUList& owner = mesh.owner();
56  const labelUList& neighbour = mesh.neighbour();
57 
58  const vectorField& C = mesh.C();
59 
60  scalarField& pLim = limiterField.primitiveFieldRef();
61 
62  forAll(pLim, face)
63  {
64  label own = owner[face];
65  label nei = neighbour[face];
66 
67  pLim[face] = Limiter::limiter
68  (
69  CDweights[face],
70  this->faceFlux_[face],
71  lPhi[own],
72  lPhi[nei],
73  gradc[own],
74  gradc[nei],
75  C[nei] - C[own]
76  );
77  }
78 
79  const typename GeometricField<Type, fvPatchField, volMesh>::Boundary&
80  bPhi = phi.boundaryField();
81 
82  surfaceScalarField::Boundary& bLim =
83  limiterField.boundaryFieldRef();
84 
85  forAll(bLim, patchi)
86  {
87  scalarField& pLim = bLim[patchi];
88 
89  if (bPhi[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 
136 
137 // * * * * * * * * * * * * Public Member Functions * * * * * * * * * * * * //
138 
139 template<class Type, class Limiter, template<class> class LimitFunc>
142 (
144 ) const
145 {
146  const fvMesh& mesh = this->mesh();
147 
148  const word limiterFieldName(type() + "Limiter(" + phi.name() + ')');
149 
150  if (this->mesh().solution().cache("limiter"))
151  {
152  if (!mesh.foundObject<surfaceScalarField>(limiterFieldName))
153  {
154  surfaceScalarField* limiterField
155  (
157  (
158  IOobject
159  (
160  limiterFieldName,
161  mesh.time().timeName(),
162  mesh,
163  IOobject::NO_READ,
164  IOobject::NO_WRITE
165  ),
166  mesh,
167  dimless
168  )
169  );
170 
171  mesh.objectRegistry::store(limiterField);
172  }
173 
174  surfaceScalarField& limiterField =
176  (
177  limiterFieldName
178  );
179 
180  calcLimiter(phi, limiterField);
181 
182  return limiterField;
183  }
184  else
185  {
186  tmp<surfaceScalarField> tlimiterField
187  (
189  (
190  limiterFieldName,
191  mesh,
192  dimless
193  )
194  );
195 
196  calcLimiter(phi, tlimiterField.ref());
197 
198  return tlimiterField;
199  }
200 }
201 
202 
203 // ************************************************************************* //
Foam::surfaceFields.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
Definition: IOobject.H:315
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Class to create NVD/TVD limited weighting-factors.
Definition: LimitedScheme.H:63
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
bool foundObject(const word &name) const
Is the named Type found?
volVectorField vectorField(fieldObject, mesh)
Generic GeometricField class.
const dimensionSet dimless
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
fvMesh & mesh
UList< label > labelUList
Definition: UList.H:65
Calculate the gradient of the given field.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
void limiter(surfaceScalarField &lambda, 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)
volScalarField scalarField(fieldObject, mesh)
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98