initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / fvm / fvmSup.C
blob67c18d3c59e9c4ab535c53199f52f23135b493cf
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 "volFields.H"
28 #include "surfaceFields.H"
29 #include "fvMatrix.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 template<class Type>
34 Foam::tmp<Foam::fvMatrix<Type> >
35 Foam::fvm::Su
37     const DimensionedField<Type, volMesh>& su,
38     GeometricField<Type, fvPatchField, volMesh>& vf
41     const fvMesh& mesh = vf.mesh();
43     tmp<fvMatrix<Type> > tfvm
44     (
45         new fvMatrix<Type>
46         (
47             vf,
48             dimVol*su.dimensions()
49         )
50     );
51     fvMatrix<Type>& fvm = tfvm();
53     fvm.source() -= mesh.V()*su.field();
55     return tfvm;
58 template<class Type>
59 Foam::tmp<Foam::fvMatrix<Type> >
60 Foam::fvm::Su
62     const tmp<DimensionedField<Type, volMesh> >& tsu,
63     GeometricField<Type, fvPatchField, volMesh>& vf
66     tmp<fvMatrix<Type> > tfvm = fvm::Su(tsu(), vf);
67     tsu.clear();
68     return tfvm;
71 template<class Type>
72 Foam::tmp<Foam::fvMatrix<Type> >
73 Foam::fvm::Su
75     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
76     GeometricField<Type, fvPatchField, volMesh>& vf
79     tmp<fvMatrix<Type> > tfvm = fvm::Su(tsu(), vf);
80     tsu.clear();
81     return tfvm;
84 template<class Type>
85 Foam::zeroField
86 Foam::fvm::Su
88     const zeroField&,
89     GeometricField<Type, fvPatchField, volMesh>& vf
92     return zeroField();
96 template<class Type>
97 Foam::tmp<Foam::fvMatrix<Type> >
98 Foam::fvm::Sp
100     const DimensionedField<scalar, volMesh>& sp,
101     GeometricField<Type, fvPatchField, volMesh>& vf
104     const fvMesh& mesh = vf.mesh();
106     tmp<fvMatrix<Type> > tfvm
107     (
108         new fvMatrix<Type>
109         (
110             vf,
111             dimVol*sp.dimensions()*vf.dimensions()
112         )
113     );
114     fvMatrix<Type>& fvm = tfvm();
116     fvm.diag() += mesh.V()*sp.field();
118     return tfvm;
121 template<class Type>
122 Foam::tmp<Foam::fvMatrix<Type> >
123 Foam::fvm::Sp
125     const tmp<DimensionedField<scalar, volMesh> >& tsp,
126     GeometricField<Type, fvPatchField, volMesh>& vf
129     tmp<fvMatrix<Type> > tfvm = fvm::Sp(tsp(), vf);
130     tsp.clear();
131     return tfvm;
134 template<class Type>
135 Foam::tmp<Foam::fvMatrix<Type> >
136 Foam::fvm::Sp
138     const tmp<volScalarField>& tsp,
139     GeometricField<Type, fvPatchField, volMesh>& vf
142     tmp<fvMatrix<Type> > tfvm = fvm::Sp(tsp(), vf);
143     tsp.clear();
144     return tfvm;
148 template<class Type>
149 Foam::tmp<Foam::fvMatrix<Type> >
150 Foam::fvm::Sp
152     const dimensionedScalar& sp,
153     GeometricField<Type, fvPatchField, volMesh>& vf
156     const fvMesh& mesh = vf.mesh();
158     tmp<fvMatrix<Type> > tfvm
159     (
160         new fvMatrix<Type>
161         (
162             vf,
163             dimVol*sp.dimensions()*vf.dimensions()
164         )
165     );
166     fvMatrix<Type>& fvm = tfvm();
168     fvm.diag() += mesh.V()*sp.value();
170     return tfvm;
173 template<class Type>
174 Foam::zeroField
175 Foam::fvm::Sp
177     const zeroField&,
178     GeometricField<Type, fvPatchField, volMesh>&
181     return zeroField();
185 template<class Type>
186 Foam::tmp<Foam::fvMatrix<Type> >
187 Foam::fvm::SuSp
189     const DimensionedField<scalar, volMesh>& susp,
190     GeometricField<Type, fvPatchField, volMesh>& vf
193     const fvMesh& mesh = vf.mesh();
195     tmp<fvMatrix<Type> > tfvm
196     (
197         new fvMatrix<Type>
198         (
199             vf,
200             dimVol*susp.dimensions()*vf.dimensions()
201         )
202     );
203     fvMatrix<Type>& fvm = tfvm();
205     fvm.diag() += mesh.V()*max(susp.field(), scalar(0));
207     fvm.source() -= mesh.V()*min(susp.field(), scalar(0))
208         *vf.internalField();
210     return tfvm;
213 template<class Type>
214 Foam::tmp<Foam::fvMatrix<Type> >
215 Foam::fvm::SuSp
217     const tmp<DimensionedField<scalar, volMesh> >& tsusp,
218     GeometricField<Type, fvPatchField, volMesh>& vf
221     tmp<fvMatrix<Type> > tfvm = fvm::SuSp(tsusp(), vf);
222     tsusp.clear();
223     return tfvm;
226 template<class Type>
227 Foam::tmp<Foam::fvMatrix<Type> >
228 Foam::fvm::SuSp
230     const tmp<volScalarField>& tsusp,
231     GeometricField<Type, fvPatchField, volMesh>& vf
234     tmp<fvMatrix<Type> > tfvm = fvm::SuSp(tsusp(), vf);
235     tsusp.clear();
236     return tfvm;
239 template<class Type>
240 Foam::zeroField
241 Foam::fvm::SuSp
243     const zeroField&,
244     GeometricField<Type, fvPatchField, volMesh>& vf
247     return zeroField();
251 // ************************************************************************* //