fvcDiv.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "fvcDiv.H"
29 #include "fvMesh.H"
30 #include "fvcSurfaceIntegrate.H"
31 #include "divScheme.H"
32 #include "convectionScheme.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace fvc
42 {
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 template<class Type>
47 tmp<GeometricField<Type, fvPatchField, volMesh>>
48 div
49 (
51 )
52 {
54  (
56  (
57  "div("+ssf.name()+')',
59  )
60  );
61 }
62 
63 
64 template<class Type>
66 div
67 (
69 )
70 {
72  tssf.clear();
73  return Div;
74 }
75 
76 
77 template<class Type>
78 tmp
79 <
81  <
83  >
84 >
85 div
86 (
88  const word& name
89 )
90 {
92  (
93  vf.mesh(), vf.mesh().divScheme(name)
94  ).ref().fvcDiv(vf);
95 }
96 
97 
98 template<class Type>
99 tmp
100 <
102  <
104  >
105 >
106 div
107 (
109  const word& name
110 )
111 {
112  typedef typename innerProduct<vector, Type>::type DivType;
114  (
115  fvc::div(tvvf(), name)
116  );
117  tvvf.clear();
118  return Div;
119 }
120 
121 template<class Type>
122 tmp
123 <
125  <
127  >
128 >
129 div
130 (
132 )
133 {
134  return fvc::div(vf, "div("+vf.name()+')');
135 }
136 
137 
138 template<class Type>
139 tmp
140 <
142  <
144  >
145 >
146 div
147 (
149 )
150 {
151  typedef typename innerProduct<vector, Type>::type DivType;
153  tvvf.clear();
154  return Div;
155 }
156 
157 
158 template<class Type>
160 div
161 (
162  const surfaceScalarField& flux,
164  const word& name
165 )
166 {
168  (
169  vf.mesh(),
170  flux,
171  vf.mesh().divScheme(name)
172  ).ref().fvcDiv(flux, vf);
173 }
174 
175 
176 template<class Type>
178 div
179 (
180  const tmp<surfaceScalarField>& tflux,
182  const word& name
183 )
184 {
186  (
187  fvc::div(tflux(), vf, name)
188  );
189  tflux.clear();
190  return Div;
191 }
192 
193 
194 template<class Type>
196 div
197 (
198  const surfaceScalarField& flux,
200  const word& name
201 )
202 {
204  (
205  fvc::div(flux, tvf(), name)
206  );
207  tvf.clear();
208  return Div;
209 }
210 
211 
212 template<class Type>
214 div
215 (
216  const tmp<surfaceScalarField>& tflux,
218  const word& name
219 )
220 {
222  (
223  fvc::div(tflux(), tvf(), name)
224  );
225  tflux.clear();
226  tvf.clear();
227  return Div;
228 }
229 
230 
231 template<class Type>
233 div
234 (
235  const surfaceScalarField& flux,
237 )
238 {
239  return fvc::div
240  (
241  flux, vf, "div("+flux.name()+','+vf.name()+')'
242  );
243 }
244 
245 
246 template<class Type>
248 div
249 (
250  const tmp<surfaceScalarField>& tflux,
252 )
253 {
255  (
256  fvc::div(tflux(), vf)
257  );
258  tflux.clear();
259  return Div;
260 }
261 
262 
263 template<class Type>
265 div
266 (
267  const surfaceScalarField& flux,
269 )
270 {
272  (
273  fvc::div(flux, tvf())
274  );
275  tvf.clear();
276  return Div;
277 }
278 
279 
280 template<class Type>
282 div
283 (
284  const tmp<surfaceScalarField>& tflux,
286 )
287 {
289  (
290  fvc::div(tflux(), tvf())
291  );
292  tflux.clear();
293  tvf.clear();
294  return Div;
295 }
296 
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 } // End namespace fvc
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 } // End namespace Foam
305 
306 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
convectionScheme.H
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:287
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::innerProduct
Definition: products.H:141
Foam::volMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:51
Foam::fv::convectionScheme::New
static tmp< convectionScheme< Type > > New(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
Return a pointer to a new convectionScheme created on freestore.
Definition: convectionScheme.C:61
fvcDiv.H
Calculate the divergence of the given field.
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
ref
rDeltaT ref()
Foam::fv::divScheme::New
static tmp< divScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new divScheme created on freestore.
Definition: divScheme.C:50
fvcSurfaceIntegrate.H
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvc::surfaceIntegrate
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcSurfaceIntegrate.C:46
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::flux
Definition: vector.H:56
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
divScheme.H