Fix tutorials: typo in tutorials/viscoelastic/viscoelasticFluidFoam/S-MDCPP/constant...
[OpenFOAM-1.6-ext.git] / src / finiteVolume / interpolation / surfaceInterpolation / limitedSchemes / limitedSurfaceInterpolationScheme / limitedSurfaceInterpolationScheme.C
blobd09ff3a43387104188f3d01997b068fc19b365ad
2 /*---------------------------------------------------------------------------*\
3   =========                 |
4   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
5    \\    /   O peration     |
6     \\  /    A nd           | Copyright held by original author
7      \\/     M anipulation  |
8 -------------------------------------------------------------------------------
9 License
10     This file is part of OpenFOAM.
12     OpenFOAM is free software; you can redistribute it and/or modify it
13     under the terms of the GNU General Public License as published by the
14     Free Software Foundation; either version 2 of the License, or (at your
15     option) any later version.
17     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
18     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
19     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
20     for more details.
22     You should have received a copy of the GNU General Public License
23     along with OpenFOAM; if not, write to the Free Software Foundation,
24     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
26 \*---------------------------------------------------------------------------*/
28 #include "limitedSurfaceInterpolationScheme.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "coupledFvPatchField.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 namespace Foam
38 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
40 template<class Type>
41 tmp<limitedSurfaceInterpolationScheme<Type> >
42 limitedSurfaceInterpolationScheme<Type>::New
44     const fvMesh& mesh,
45     Istream& schemeData
48     if (surfaceInterpolation::debug)
49     {
50         Info<< "limitedSurfaceInterpolationScheme<Type>::"
51                "New(const fvMesh&, Istream&)"
52                " : constructing limitedSurfaceInterpolationScheme<Type>"
53             << endl;
54     }
56     if (schemeData.eof())
57     {
58         FatalIOErrorIn
59         (
60             "limitedSurfaceInterpolationScheme<Type>::"
61             "New(const fvMesh&, Istream&)",
62             schemeData
63         )   << "Discretisation scheme not specified"
64             << endl << endl
65             << "Valid schemes are :" << endl
66             << MeshConstructorTablePtr_->toc()
67             << exit(FatalIOError);
68     }
70     word schemeName(schemeData);
72     typename MeshConstructorTable::iterator constructorIter =
73         MeshConstructorTablePtr_->find(schemeName);
75     if (constructorIter == MeshConstructorTablePtr_->end())
76     {
77         FatalIOErrorIn
78         (
79             "limitedSurfaceInterpolationScheme<Type>::"
80             "New(const fvMesh&, Istream&)",
81             schemeData
82         )   << "Unknown discretisation scheme " << schemeName
83             << endl << endl
84             << "Valid schemes are :" << endl
85             << MeshConstructorTablePtr_->toc()
86             << exit(FatalIOError);
87     }
89     return constructorIter()(mesh, schemeData);
93 // Return weighting factors for scheme given by name in dictionary
94 template<class Type>
95 tmp<limitedSurfaceInterpolationScheme<Type> >
96 limitedSurfaceInterpolationScheme<Type>::New
98     const fvMesh& mesh,
99     const surfaceScalarField& faceFlux,
100     Istream& schemeData
103     if (surfaceInterpolation::debug)
104     {
105         Info<< "limitedSurfaceInterpolationScheme<Type>::New"
106                "(const fvMesh&, const surfaceScalarField&, Istream&) : "
107                "constructing limitedSurfaceInterpolationScheme<Type>"
108             << endl;
109     }
111     if (schemeData.eof())
112     {
113         FatalIOErrorIn
114         (
115             "limitedSurfaceInterpolationScheme<Type>::New"
116             "(const fvMesh&, const surfaceScalarField&, Istream&)",
117             schemeData
118         )   << "Discretisation scheme not specified"
119             << endl << endl
120             << "Valid schemes are :" << endl
121             << MeshConstructorTablePtr_->toc()
122             << exit(FatalIOError);
123     }
125     word schemeName(schemeData);
127     typename MeshFluxConstructorTable::iterator constructorIter =
128         MeshFluxConstructorTablePtr_->find(schemeName);
130     if (constructorIter == MeshFluxConstructorTablePtr_->end())
131     {
132         FatalIOErrorIn
133         (
134             "limitedSurfaceInterpolationScheme<Type>::New"
135             "(const fvMesh&, const surfaceScalarField&, Istream&)",
136             schemeData
137         )   << "Unknown discretisation scheme " << schemeName
138             << endl << endl
139             << "Valid schemes are :" << endl
140             << MeshFluxConstructorTablePtr_->toc()
141             << exit(FatalIOError);
142     }
144     return constructorIter()(mesh, faceFlux, schemeData);
148 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
150 template<class Type>
151 limitedSurfaceInterpolationScheme<Type>::~limitedSurfaceInterpolationScheme()
155 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
157 template<class Type>
158 tmp<surfaceScalarField> limitedSurfaceInterpolationScheme<Type>::weights
160     const GeometricField<Type, fvPatchField, volMesh>& phi,
161     const surfaceScalarField& CDweights,
162     tmp<surfaceScalarField> tLimiter
163 ) const
165     // Note that here the weights field is initialised as the limiter
166     // from which the weight is calculated using the limiter value
167     surfaceScalarField& Weights = tLimiter();
169     scalarField& pWeights = Weights.internalField();
171     forAll(pWeights, face)
172     {
173         pWeights[face] =
174             pWeights[face]*CDweights[face]
175           + (1.0 - pWeights[face])*pos(faceFlux_[face]);
176     }
178     surfaceScalarField::GeometricBoundaryField& bWeights =
179         Weights.boundaryField();
181     forAll(bWeights, patchI)
182     {
183         scalarField& pWeights = bWeights[patchI];
185         const scalarField& pCDweights = CDweights.boundaryField()[patchI];
186         const scalarField& pFaceFlux = faceFlux_.boundaryField()[patchI];
188         forAll(pWeights, face)
189         {
190             pWeights[face] =
191                 pWeights[face]*pCDweights[face]
192               + (1.0 - pWeights[face])*pos(pFaceFlux[face]);
193         }
194     }
196     return tLimiter;
199 template<class Type>
200 tmp<surfaceScalarField> limitedSurfaceInterpolationScheme<Type>::weights
202     const GeometricField<Type, fvPatchField, volMesh>& phi
203 ) const
205     return this->weights
206     (
207         phi,
208         this->mesh().surfaceInterpolation::weights(),
209         this->limiter(phi)
210     );
213 template<class Type>
214 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
215 limitedSurfaceInterpolationScheme<Type>::flux
217     const GeometricField<Type, fvPatchField, volMesh>& phi
218 ) const
220     return faceFlux_*interpolate(phi);
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 } // End namespace Foam
228 // ************************************************************************* //