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 VolField<Type>& phi,
37  surfaceScalarField& limiterField
38 ) const
39 {
40  const fvMesh& mesh = this->mesh();
41 
42  tmp<VolField<typename Limiter::phiType>> tlPhi = LimitFunc<Type>()(phi);
43  const VolField<typename Limiter::phiType>& lPhi = tlPhi();
44 
45  tmp<VolField<typename Limiter::gradPhiType>> tgradc(fvc::grad(lPhi));
46  const VolField<typename Limiter::gradPhiType>& gradc = tgradc();
47 
48  const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights();
49 
50  const labelUList& owner = mesh.owner();
51  const labelUList& neighbour = mesh.neighbour();
52 
53  const vectorField& C = mesh.C();
54 
55  scalarField& pLim = limiterField.primitiveFieldRef();
56 
57  forAll(pLim, face)
58  {
59  label own = owner[face];
60  label nei = neighbour[face];
61 
62  pLim[face] = Limiter::limiter
63  (
64  CDweights[face],
65  this->faceFlux_[face],
66  lPhi[own],
67  lPhi[nei],
68  gradc[own],
69  gradc[nei],
70  C[nei] - C[own]
71  );
72  }
73 
74  const typename VolField<Type>::Boundary&
75  bPhi = phi.boundaryField();
76 
77  surfaceScalarField::Boundary& bLim =
78  limiterField.boundaryFieldRef();
79 
80  forAll(bLim, patchi)
81  {
82  scalarField& pLim = bLim[patchi];
83 
84  if (bPhi[patchi].coupled())
85  {
86  const scalarField& pCDweights = CDweights.boundaryField()[patchi];
87  const scalarField& pFaceFlux =
88  this->faceFlux_.boundaryField()[patchi];
89 
90  const Field<typename Limiter::phiType> plPhiP
91  (
92  lPhi.boundaryField()[patchi].patchInternalField()
93  );
94  const Field<typename Limiter::phiType> plPhiN
95  (
96  lPhi.boundaryField()[patchi].patchNeighbourField()
97  );
98  const Field<typename Limiter::gradPhiType> pGradcP
99  (
100  gradc.boundaryField()[patchi].patchInternalField()
101  );
102  const Field<typename Limiter::gradPhiType> pGradcN
103  (
104  gradc.boundaryField()[patchi].patchNeighbourField()
105  );
106 
107  // Build the d-vectors
108  vectorField pd(CDweights.boundaryField()[patchi].patch().delta());
109 
110  forAll(pLim, face)
111  {
112  pLim[face] = Limiter::limiter
113  (
114  pCDweights[face],
115  pFaceFlux[face],
116  plPhiP[face],
117  plPhiN[face],
118  pGradcP[face],
119  pGradcN[face],
120  pd[face]
121  );
122  }
123  }
124  else
125  {
126  pLim = 1.0;
127  }
128  }
129 }
130 
131 
132 // * * * * * * * * * * * * Public Member Functions * * * * * * * * * * * * //
133 
134 template<class Type, class Limiter, template<class> class LimitFunc>
137 (
138  const VolField<Type>& phi
139 ) const
140 {
141  const fvMesh& mesh = this->mesh();
142 
143  const word limiterFieldName(type() + "Limiter(" + phi.name() + ')');
144 
145  if (this->mesh().solution().cache("limiter"))
146  {
147  if (!mesh.foundObject<surfaceScalarField>(limiterFieldName))
148  {
149  surfaceScalarField* limiterField
150  (
152  (
153  IOobject
154  (
155  limiterFieldName,
156  mesh.time().name(),
157  mesh,
158  IOobject::NO_READ,
159  IOobject::NO_WRITE
160  ),
161  mesh,
162  dimless
163  )
164  );
165 
166  mesh.objectRegistry::store(limiterField);
167  }
168 
169  surfaceScalarField& limiterField =
171  (
172  limiterFieldName
173  );
174 
175  calcLimiter(phi, limiterField);
176 
177  return limiterField;
178  }
179  else
180  {
181  tmp<surfaceScalarField> tlimiterField
182  (
184  (
185  limiterFieldName,
186  mesh,
187  dimless
188  )
189  );
190 
191  calcLimiter(phi, tlimiterField.ref());
192 
193  return tlimiterField;
194  }
195 }
196 
197 
198 // ************************************************************************* //
static const Foam::dimensionedScalar C("C", Foam::dimTemperature, 234.5)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
Class to create NVD/TVD limited weighting-factors.
Definition: LimitedScheme.H:67
virtual tmp< surfaceScalarField > limiter(const VolField< Type > &) const
Return the interpolation weighting factors.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
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 in registry.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
Calculate the gradient of the given field.
volScalarField scalarField(fieldObject, mesh)
volVectorField vectorField(fieldObject, mesh)
label patchi
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)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
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:59
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
UList< label > labelUList
Definition: UList.H:65
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Foam::surfaceFields.