limitedSurfaceInterpolationScheme.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 | Copyright (C) 2011-2015 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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
37 
38 template<class Type>
39 tmp<limitedSurfaceInterpolationScheme<Type> >
41 (
42  const fvMesh& mesh,
43  Istream& schemeData
44 )
45 {
46  if (surfaceInterpolation::debug)
47  {
48  Info<< "limitedSurfaceInterpolationScheme<Type>::"
49  "New(const fvMesh&, Istream&)"
50  " : constructing limitedSurfaceInterpolationScheme<Type>"
51  << endl;
52  }
53 
54  if (schemeData.eof())
55  {
57  (
58  schemeData
59  ) << "Discretisation scheme not specified"
60  << endl << endl
61  << "Valid schemes are :" << endl
62  << MeshConstructorTablePtr_->sortedToc()
63  << exit(FatalIOError);
64  }
65 
66  const word schemeName(schemeData);
67 
68  typename MeshConstructorTable::iterator constructorIter =
69  MeshConstructorTablePtr_->find(schemeName);
70 
71  if (constructorIter == MeshConstructorTablePtr_->end())
72  {
74  (
75  schemeData
76  ) << "Unknown discretisation scheme "
77  << schemeName << nl << nl
78  << "Valid schemes are :" << endl
79  << MeshConstructorTablePtr_->sortedToc()
80  << exit(FatalIOError);
81  }
82 
83  return constructorIter()(mesh, schemeData);
84 }
85 
86 
87 // Return weighting factors for scheme given by name in dictionary
88 template<class Type>
91 (
92  const fvMesh& mesh,
93  const surfaceScalarField& faceFlux,
94  Istream& schemeData
95 )
96 {
97  if (surfaceInterpolation::debug)
98  {
99  Info<< "limitedSurfaceInterpolationScheme<Type>::New"
100  "(const fvMesh&, const surfaceScalarField&, Istream&) : "
101  "constructing limitedSurfaceInterpolationScheme<Type>"
102  << endl;
103  }
104 
105  if (schemeData.eof())
106  {
108  (
109  schemeData
110  ) << "Discretisation scheme not specified"
111  << endl << endl
112  << "Valid schemes are :" << endl
113  << MeshConstructorTablePtr_->sortedToc()
114  << exit(FatalIOError);
115  }
116 
117  const word schemeName(schemeData);
118 
119  typename MeshFluxConstructorTable::iterator constructorIter =
120  MeshFluxConstructorTablePtr_->find(schemeName);
121 
122  if (constructorIter == MeshFluxConstructorTablePtr_->end())
123  {
125  (
126  schemeData
127  ) << "Unknown discretisation scheme "
128  << schemeName << nl << nl
129  << "Valid schemes are :" << endl
130  << MeshFluxConstructorTablePtr_->sortedToc()
131  << exit(FatalIOError);
132  }
133 
134  return constructorIter()(mesh, faceFlux, schemeData);
135 }
136 
137 
138 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
139 
140 template<class Type>
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
147 template<class Type>
149 (
151  const surfaceScalarField& CDweights,
152  tmp<surfaceScalarField> tLimiter
153 ) const
154 {
155  // Note that here the weights field is initialised as the limiter
156  // from which the weight is calculated using the limiter value
157  surfaceScalarField& Weights = tLimiter();
158 
159  scalarField& pWeights = Weights.internalField();
160 
161  forAll(pWeights, face)
162  {
163  pWeights[face] =
164  pWeights[face]*CDweights[face]
165  + (1.0 - pWeights[face])*pos(faceFlux_[face]);
166  }
167 
169  Weights.boundaryField();
170 
171  forAll(bWeights, patchI)
172  {
173  scalarField& pWeights = bWeights[patchI];
174 
175  const scalarField& pCDweights = CDweights.boundaryField()[patchI];
176  const scalarField& pFaceFlux = faceFlux_.boundaryField()[patchI];
177 
178  forAll(pWeights, face)
179  {
180  pWeights[face] =
181  pWeights[face]*pCDweights[face]
182  + (1.0 - pWeights[face])*pos(pFaceFlux[face]);
183  }
184  }
185 
186  return tLimiter;
187 }
188 
189 template<class Type>
191 (
193 ) const
194 {
195  return this->weights
196  (
197  phi,
199  this->limiter(phi)
200  );
201 }
202 
203 template<class Type>
206 (
208 ) const
209 {
210  return faceFlux_*this->interpolate(phi);
211 }
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 } // End namespace Foam
217 
218 // ************************************************************************* //
volFields.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::IOstream::eof
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:339
Foam::limitedSurfaceInterpolationScheme::~limitedSurfaceInterpolationScheme
virtual ~limitedSurfaceInterpolationScheme()
Destructor.
Definition: limitedSurfaceInterpolationScheme.C:141
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::limiter
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 scalar psiMax, const scalar psiMin)
Definition: MULESTemplates.C:159
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::surfaceInterpolation::weights
const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
Definition: surfaceInterpolation.C:77
Foam::limitedSurfaceInterpolationScheme::New
static tmp< limitedSurfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Definition: limitedSurfaceInterpolationScheme.C:41
limitedSurfaceInterpolationScheme.H
Foam::limitedSurfaceInterpolationScheme::flux
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Definition: limitedSurfaceInterpolationScheme.C:206
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
coupledFvPatchField.H
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::GeometricField::GeometricBoundaryField
Definition: GeometricField.H:105
Foam::limitedSurfaceInterpolationScheme::weights
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &, const surfaceScalarField &CDweights, tmp< surfaceScalarField > tLimiter) const
Return the interpolation weighting factors for the given field,.
Definition: limitedSurfaceInterpolationScheme.C:149
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190