fvcDiv.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvcDiv.H"
27 #include "fvMesh.H"
28 #include "fvcSurfaceIntegrate.H"
29 #include "divScheme.H"
30 #include "convectionScheme.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fvc
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<class Type>
45 tmp<GeometricField<Type, fvPatchField, volMesh>>
46 div
47 (
49 )
50 {
52  (
54  (
55  "div("+ssf.name()+')',
57  )
58  );
59 }
60 
61 
62 template<class Type>
64 div
65 (
67 )
68 {
70  tssf.clear();
71  return Div;
72 }
73 
74 
75 template<class Type>
76 tmp
77 <
79  <
80  typename innerProduct<vector, Type>::type, fvPatchField, volMesh
81  >
82 >
83 div
84 (
86  const word& name
87 )
88 {
90  (
91  vf.mesh(), vf.mesh().divScheme(name)
92  ).ref().fvcDiv(vf);
93 }
94 
95 
96 template<class Type>
97 tmp
98 <
100  <
101  typename innerProduct<vector, Type>::type, fvPatchField, volMesh
102  >
103 >
104 div
105 (
107  const word& name
108 )
109 {
110  typedef typename innerProduct<vector, Type>::type DivType;
112  (
113  fvc::div(tvvf(), name)
114  );
115  tvvf.clear();
116  return Div;
117 }
118 
119 template<class Type>
120 tmp
121 <
123  <
124  typename innerProduct<vector, Type>::type, fvPatchField, volMesh
125  >
126 >
127 div
128 (
130 )
131 {
132  return fvc::div(vf, "div("+vf.name()+')');
133 }
134 
135 
136 template<class Type>
137 tmp
138 <
140  <
141  typename innerProduct<vector, Type>::type, fvPatchField, volMesh
142  >
143 >
144 div
145 (
147 )
148 {
149  typedef typename innerProduct<vector, Type>::type DivType;
151  tvvf.clear();
152  return Div;
153 }
154 
155 
156 template<class Type>
158 div
159 (
160  const surfaceScalarField& flux,
162  const word& name
163 )
164 {
166  (
167  vf.mesh(),
168  flux,
169  vf.mesh().divScheme(name)
170  ).ref().fvcDiv(flux, vf);
171 }
172 
173 
174 template<class Type>
176 div
177 (
178  const tmp<surfaceScalarField>& tflux,
180  const word& name
181 )
182 {
184  (
185  fvc::div(tflux(), vf, name)
186  );
187  tflux.clear();
188  return Div;
189 }
190 
191 
192 template<class Type>
194 div
195 (
196  const surfaceScalarField& flux,
198  const word& name
199 )
200 {
202  (
203  fvc::div(flux, tvf(), name)
204  );
205  tvf.clear();
206  return Div;
207 }
208 
209 
210 template<class Type>
212 div
213 (
214  const tmp<surfaceScalarField>& tflux,
216  const word& name
217 )
218 {
220  (
221  fvc::div(tflux(), tvf(), name)
222  );
223  tflux.clear();
224  tvf.clear();
225  return Div;
226 }
227 
228 
229 template<class Type>
231 div
232 (
233  const surfaceScalarField& flux,
235 )
236 {
237  return fvc::div
238  (
239  flux, vf, "div("+flux.name()+','+vf.name()+')'
240  );
241 }
242 
243 
244 template<class Type>
246 div
247 (
248  const tmp<surfaceScalarField>& tflux,
250 )
251 {
253  (
254  fvc::div(tflux(), vf)
255  );
256  tflux.clear();
257  return Div;
258 }
259 
260 
261 template<class Type>
263 div
264 (
265  const surfaceScalarField& flux,
267 )
268 {
270  (
271  fvc::div(flux, tvf())
272  );
273  tvf.clear();
274  return Div;
275 }
276 
277 
278 template<class Type>
280 div
281 (
282  const tmp<surfaceScalarField>& tflux,
284 )
285 {
287  (
288  fvc::div(tflux(), tvf())
289  );
290  tflux.clear();
291  tvf.clear();
292  return Div;
293 }
294 
295 
296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297 
298 } // End namespace fvc
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 } // End namespace Foam
303 
304 // ************************************************************************* //
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
const word & name() const
Return name.
Definition: IOobject.H:297
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Generic GeometricField class.
static tmp< divScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new divScheme created on freestore.
Definition: divScheme.C:47
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:53
A class for handling words, derived from string.
Definition: word.H:59
Calculate the divergence of the given field.
const Mesh & mesh() const
Return mesh.
static tmp< convectionScheme< Type > > New(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
Return a pointer to a new convectionScheme created on freestore.
rDeltaT ref()
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
Namespace for OpenFOAM.