1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 Foam::tmp<Foam::fvMatrix<Type> >
37 const DimensionedField<Type, volMesh>& su,
38 GeometricField<Type, fvPatchField, volMesh>& vf
41 const fvMesh& mesh = vf.mesh();
43 tmp<fvMatrix<Type> > tfvm
48 dimVol*su.dimensions()
51 fvMatrix<Type>& fvm = tfvm();
53 fvm.source() -= mesh.V()*su.field();
59 Foam::tmp<Foam::fvMatrix<Type> >
62 const tmp<DimensionedField<Type, volMesh> >& tsu,
63 GeometricField<Type, fvPatchField, volMesh>& vf
66 tmp<fvMatrix<Type> > tfvm = fvm::Su(tsu(), vf);
72 Foam::tmp<Foam::fvMatrix<Type> >
75 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
76 GeometricField<Type, fvPatchField, volMesh>& vf
79 tmp<fvMatrix<Type> > tfvm = fvm::Su(tsu(), vf);
89 GeometricField<Type, fvPatchField, volMesh>& vf
97 Foam::tmp<Foam::fvMatrix<Type> >
100 const DimensionedField<scalar, volMesh>& sp,
101 GeometricField<Type, fvPatchField, volMesh>& vf
104 const fvMesh& mesh = vf.mesh();
106 tmp<fvMatrix<Type> > tfvm
111 dimVol*sp.dimensions()*vf.dimensions()
114 fvMatrix<Type>& fvm = tfvm();
116 fvm.diag() += mesh.V()*sp.field();
122 Foam::tmp<Foam::fvMatrix<Type> >
125 const tmp<DimensionedField<scalar, volMesh> >& tsp,
126 GeometricField<Type, fvPatchField, volMesh>& vf
129 tmp<fvMatrix<Type> > tfvm = fvm::Sp(tsp(), vf);
135 Foam::tmp<Foam::fvMatrix<Type> >
138 const tmp<volScalarField>& tsp,
139 GeometricField<Type, fvPatchField, volMesh>& vf
142 tmp<fvMatrix<Type> > tfvm = fvm::Sp(tsp(), vf);
149 Foam::tmp<Foam::fvMatrix<Type> >
152 const dimensionedScalar& sp,
153 GeometricField<Type, fvPatchField, volMesh>& vf
156 const fvMesh& mesh = vf.mesh();
158 tmp<fvMatrix<Type> > tfvm
163 dimVol*sp.dimensions()*vf.dimensions()
166 fvMatrix<Type>& fvm = tfvm();
168 fvm.diag() += mesh.V()*sp.value();
178 GeometricField<Type, fvPatchField, volMesh>&
186 Foam::tmp<Foam::fvMatrix<Type> >
189 const DimensionedField<scalar, volMesh>& susp,
190 GeometricField<Type, fvPatchField, volMesh>& vf
193 const fvMesh& mesh = vf.mesh();
195 tmp<fvMatrix<Type> > tfvm
200 dimVol*susp.dimensions()*vf.dimensions()
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))
214 Foam::tmp<Foam::fvMatrix<Type> >
217 const tmp<DimensionedField<scalar, volMesh> >& tsusp,
218 GeometricField<Type, fvPatchField, volMesh>& vf
221 tmp<fvMatrix<Type> > tfvm = fvm::SuSp(tsusp(), vf);
227 Foam::tmp<Foam::fvMatrix<Type> >
230 const tmp<volScalarField>& tsusp,
231 GeometricField<Type, fvPatchField, volMesh>& vf
234 tmp<fvMatrix<Type> > tfvm = fvm::SuSp(tsusp(), vf);
244 GeometricField<Type, fvPatchField, volMesh>& vf
251 // ************************************************************************* //