initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / interpolation / surfaceInterpolation / limitedSchemes / limitedSurfaceInterpolationScheme / limitedSurfaceInterpolationScheme.C
blob77b18ab8066be7911b90404d7b21a70a9ec825f1
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "limitedSurfaceInterpolationScheme.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "coupledFvPatchField.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
39 template<class Type>
40 tmp<limitedSurfaceInterpolationScheme<Type> >
41 limitedSurfaceInterpolationScheme<Type>::New
43     const fvMesh& mesh,
44     Istream& schemeData
47     if (surfaceInterpolation::debug)
48     {
49         Info<< "limitedSurfaceInterpolationScheme<Type>::"
50                "New(const fvMesh&, Istream&)"
51                " : constructing limitedSurfaceInterpolationScheme<Type>"
52             << endl;
53     }
55     if (schemeData.eof())
56     {
57         FatalIOErrorIn
58         (
59             "limitedSurfaceInterpolationScheme<Type>::"
60             "New(const fvMesh&, Istream&)",
61             schemeData
62         )   << "Discretisation scheme not specified"
63             << endl << endl
64             << "Valid schemes are :" << endl
65             << MeshConstructorTablePtr_->toc()
66             << exit(FatalIOError);
67     }
69     word schemeName(schemeData);
71     typename MeshConstructorTable::iterator constructorIter =
72         MeshConstructorTablePtr_->find(schemeName);
74     if (constructorIter == MeshConstructorTablePtr_->end())
75     {
76         FatalIOErrorIn
77         (
78             "limitedSurfaceInterpolationScheme<Type>::"
79             "New(const fvMesh&, Istream&)",
80             schemeData
81         )   << "Unknown discretisation scheme " << schemeName
82             << endl << endl
83             << "Valid schemes are :" << endl
84             << MeshConstructorTablePtr_->toc()
85             << exit(FatalIOError);
86     }
88     return constructorIter()(mesh, schemeData);
92 // Return weighting factors for scheme given by name in dictionary
93 template<class Type>
94 tmp<limitedSurfaceInterpolationScheme<Type> >
95 limitedSurfaceInterpolationScheme<Type>::New
97     const fvMesh& mesh,
98     const surfaceScalarField& faceFlux,
99     Istream& schemeData
102     if (surfaceInterpolation::debug)
103     {
104         Info<< "limitedSurfaceInterpolationScheme<Type>::New"
105                "(const fvMesh&, const surfaceScalarField&, Istream&) : "
106                "constructing limitedSurfaceInterpolationScheme<Type>"
107             << endl;
108     }
110     if (schemeData.eof())
111     {
112         FatalIOErrorIn
113         (
114             "limitedSurfaceInterpolationScheme<Type>::New"
115             "(const fvMesh&, const surfaceScalarField&, Istream&)",
116             schemeData
117         )   << "Discretisation scheme not specified"
118             << endl << endl
119             << "Valid schemes are :" << endl
120             << MeshConstructorTablePtr_->toc()
121             << exit(FatalIOError);
122     }
124     word schemeName(schemeData);
126     typename MeshFluxConstructorTable::iterator constructorIter =
127         MeshFluxConstructorTablePtr_->find(schemeName);
129     if (constructorIter == MeshFluxConstructorTablePtr_->end())
130     {
131         FatalIOErrorIn
132         (
133             "limitedSurfaceInterpolationScheme<Type>::New"
134             "(const fvMesh&, const surfaceScalarField&, Istream&)",
135             schemeData
136         )   << "Unknown discretisation scheme " << schemeName
137             << endl << endl
138             << "Valid schemes are :" << endl
139             << MeshFluxConstructorTablePtr_->toc()
140             << exit(FatalIOError);
141     }
143     return constructorIter()(mesh, faceFlux, schemeData);
147 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
149 template<class Type>
150 limitedSurfaceInterpolationScheme<Type>::~limitedSurfaceInterpolationScheme()
154 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
156 template<class Type>
157 tmp<surfaceScalarField> limitedSurfaceInterpolationScheme<Type>::weights
159     const GeometricField<Type, fvPatchField, volMesh>& phi,
160     const surfaceScalarField& CDweights,
161     tmp<surfaceScalarField> tLimiter
162 ) const
164     // Note that here the weights field is initialised as the limiter
165     // from which the weight is calculated using the limiter value
166     surfaceScalarField& Weights = tLimiter();
168     scalarField& pWeights = Weights.internalField();
170     forAll(pWeights, face)
171     {
172         pWeights[face] =
173             pWeights[face]*CDweights[face]
174           + (1.0 - pWeights[face])*pos(faceFlux_[face]);
175     }
177     surfaceScalarField::GeometricBoundaryField& bWeights =
178         Weights.boundaryField();
180     forAll(bWeights, patchI)
181     {
182         scalarField& pWeights = bWeights[patchI];
184         const scalarField& pCDweights = CDweights.boundaryField()[patchI];
185         const scalarField& pFaceFlux = faceFlux_.boundaryField()[patchI];
187         forAll(pWeights, face)
188         {
189             pWeights[face] =
190                 pWeights[face]*pCDweights[face]
191               + (1.0 - pWeights[face])*pos(pFaceFlux[face]);
192         }
193     }
195     return tLimiter;
198 template<class Type>
199 tmp<surfaceScalarField> limitedSurfaceInterpolationScheme<Type>::weights
201     const GeometricField<Type, fvPatchField, volMesh>& phi
202 ) const
204     return this->weights
205     (
206         phi,
207         this->mesh().surfaceInterpolation::weights(),
208         this->limiter(phi)
209     );
212 template<class Type>
213 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
214 limitedSurfaceInterpolationScheme<Type>::flux
216     const GeometricField<Type, fvPatchField, volMesh>& phi
217 ) const
219     return faceFlux_*interpolate(phi);
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 } // End namespace Foam
227 // ************************************************************************* //