limitedSurfaceInterpolationScheme.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-2018 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 
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "coupledFvPatchField.H"
30 
31 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
32 
33 template<class Type>
36 (
37  const fvMesh& mesh,
38  Istream& schemeData
39 )
40 {
41  if (surfaceInterpolation::debug)
42  {
44  << "Constructing limitedSurfaceInterpolationScheme<Type>" << endl;
45  }
46 
47  if (schemeData.eof())
48  {
50  (
51  schemeData
52  ) << "Discretisation scheme not specified"
53  << endl << endl
54  << "Valid schemes are :" << endl
55  << MeshConstructorTablePtr_->sortedToc()
56  << exit(FatalIOError);
57  }
58 
59  const word schemeName(schemeData);
60 
61  typename MeshConstructorTable::iterator constructorIter =
62  MeshConstructorTablePtr_->find(schemeName);
63 
64  if (constructorIter == MeshConstructorTablePtr_->end())
65  {
67  (
68  schemeData
69  ) << "Unknown discretisation scheme "
70  << schemeName << nl << nl
71  << "Valid schemes are :" << endl
72  << MeshConstructorTablePtr_->sortedToc()
73  << exit(FatalIOError);
74  }
75 
76  return constructorIter()(mesh, schemeData);
77 }
78 
79 
80 template<class Type>
83 (
84  const fvMesh& mesh,
85  const surfaceScalarField& faceFlux,
86  Istream& schemeData
87 )
88 {
89  if (surfaceInterpolation::debug)
90  {
92  << "Constructing limitedSurfaceInterpolationScheme<Type>"
93  << endl;
94  }
95 
96  if (schemeData.eof())
97  {
99  (
100  schemeData
101  ) << "Discretisation scheme not specified"
102  << endl << endl
103  << "Valid schemes are :" << endl
104  << MeshConstructorTablePtr_->sortedToc()
105  << exit(FatalIOError);
106  }
107 
108  const word schemeName(schemeData);
109 
110  typename MeshFluxConstructorTable::iterator constructorIter =
111  MeshFluxConstructorTablePtr_->find(schemeName);
112 
113  if (constructorIter == MeshFluxConstructorTablePtr_->end())
114  {
116  (
117  schemeData
118  ) << "Unknown discretisation scheme "
119  << schemeName << nl << nl
120  << "Valid schemes are :" << endl
121  << MeshFluxConstructorTablePtr_->sortedToc()
122  << exit(FatalIOError);
123  }
124 
125  return constructorIter()(mesh, faceFlux, schemeData);
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
130 
131 template<class Type>
134 {}
135 
136 
137 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 
139 template<class Type>
142 (
144  const surfaceScalarField& CDweights,
145  tmp<surfaceScalarField> tLimiter
146 ) const
147 {
148  // Note that here the weights field is initialised as the limiter
149  // from which the weight is calculated using the limiter value
150  surfaceScalarField& Weights = tLimiter.ref();
151 
152  scalarField& pWeights = Weights.primitiveFieldRef();
153 
154  forAll(pWeights, face)
155  {
156  pWeights[face] =
157  pWeights[face]*CDweights[face]
158  + (1.0 - pWeights[face])*pos0(faceFlux_[face]);
159  }
160 
161  surfaceScalarField::Boundary& bWeights =
162  Weights.boundaryFieldRef();
163 
164  forAll(bWeights, patchi)
165  {
166  scalarField& pWeights = bWeights[patchi];
167 
168  const scalarField& pCDweights = CDweights.boundaryField()[patchi];
169  const scalarField& pFaceFlux = faceFlux_.boundaryField()[patchi];
170 
171  forAll(pWeights, face)
172  {
173  pWeights[face] =
174  pWeights[face]*pCDweights[face]
175  + (1.0 - pWeights[face])*pos0(pFaceFlux[face]);
176  }
177  }
178 
179  return tLimiter;
180 }
181 
182 template<class Type>
185 (
187 ) const
188 {
189  return this->weights
190  (
191  phi,
192  this->mesh().surfaceInterpolation::weights(),
193  this->limiter(phi)
194  );
195 }
196 
197 template<class Type>
200 (
202 ) const
203 {
204  return faceFlux_*this->interpolate(phi);
205 }
206 
207 
208 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
static tmp< limitedSurfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
dynamicFvMesh & mesh
void limiter(scalarField &allLambda, 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)
A class for handling words, derived from string.
Definition: word.H:59
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:339
dimensionedScalar pos0(const dimensionedScalar &ds)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &, const surfaceScalarField &CDweights, tmp< surfaceScalarField > tLimiter) const
Return the interpolation weighting factors for the given field,.
static const char nl
Definition: Ostream.H:265
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
label patchi
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:53
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.