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 \*---------------------------------------------------------------------------*/
29 #include "fvcSurfaceIntegrate.H"
30 #include "divScheme.H"
31 #include "convectionScheme.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 tmp<GeometricField<Type, fvPatchField, volMesh> >
49 const GeometricField<Type, fvsPatchField, surfaceMesh>& ssf
52 return tmp<GeometricField<Type, fvPatchField, volMesh> >
54 new GeometricField<Type, fvPatchField, volMesh>
56 "div("+ssf.name()+')',
57 fvc::surfaceIntegrate(ssf)
64 tmp<GeometricField<Type, fvPatchField, volMesh> >
67 const tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >& tssf
70 tmp<GeometricField<Type, fvPatchField, volMesh> > Div(fvc::div(tssf()));
81 typename innerProduct<vector, Type>::type, fvPatchField, volMesh
86 const GeometricField<Type, fvPatchField, volMesh>& vf,
90 return fv::divScheme<Type>::New
92 vf.mesh(), vf.mesh().divScheme(name)
102 typename innerProduct<vector, Type>::type, fvPatchField, volMesh
107 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvvf,
111 typedef typename innerProduct<vector, Type>::type DivType;
112 tmp<GeometricField<DivType, fvPatchField, volMesh> > Div
114 fvc::div(tvvf(), name)
125 typename innerProduct<vector, Type>::type, fvPatchField, volMesh
130 const GeometricField<Type, fvPatchField, volMesh>& vf
133 return fvc::div(vf, "div("+vf.name()+')');
142 typename innerProduct<vector, Type>::type, fvPatchField, volMesh
147 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvvf
150 typedef typename innerProduct<vector, Type>::type DivType;
151 tmp<GeometricField<DivType, fvPatchField, volMesh> > Div(fvc::div(tvvf()));
158 tmp<GeometricField<Type, fvPatchField, volMesh> >
161 const surfaceScalarField& flux,
162 const GeometricField<Type, fvPatchField, volMesh>& vf,
166 return fv::convectionScheme<Type>::New
170 vf.mesh().divScheme(name)
171 )().fvcDiv(flux, vf);
176 tmp<GeometricField<Type, fvPatchField, volMesh> >
179 const tmp<surfaceScalarField>& tflux,
180 const GeometricField<Type, fvPatchField, volMesh>& vf,
184 tmp<GeometricField<Type, fvPatchField, volMesh> > Div
186 fvc::div(tflux(), vf, name)
194 tmp<GeometricField<Type, fvPatchField, volMesh> >
197 const surfaceScalarField& flux,
198 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
202 tmp<GeometricField<Type, fvPatchField, volMesh> > Div
204 fvc::div(flux, tvf(), name)
212 tmp<GeometricField<Type, fvPatchField, volMesh> >
215 const tmp<surfaceScalarField>& tflux,
216 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
220 tmp<GeometricField<Type, fvPatchField, volMesh> > Div
222 fvc::div(tflux(), tvf(), name)
231 tmp<GeometricField<Type, fvPatchField, volMesh> >
234 const surfaceScalarField& flux,
235 const GeometricField<Type, fvPatchField, volMesh>& vf
240 flux, vf, "div("+flux.name()+','+vf.name()+')'
246 tmp<GeometricField<Type, fvPatchField, volMesh> >
249 const tmp<surfaceScalarField>& tflux,
250 const GeometricField<Type, fvPatchField, volMesh>& vf
253 tmp<GeometricField<Type, fvPatchField, volMesh> > Div
255 fvc::div(tflux(), vf)
263 tmp<GeometricField<Type, fvPatchField, volMesh> >
266 const surfaceScalarField& flux,
267 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf
270 tmp<GeometricField<Type, fvPatchField, volMesh> > Div
272 fvc::div(flux, tvf())
280 tmp<GeometricField<Type, fvPatchField, volMesh> >
283 const tmp<surfaceScalarField>& tflux,
284 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf
287 tmp<GeometricField<Type, fvPatchField, volMesh> > Div
289 fvc::div(tflux(), tvf())
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 } // End namespace fvc
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 } // End namespace Foam
305 // ************************************************************************* //