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-2022 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  (
53  "div("+ssf.name()+')',
55  );
56 }
57 
58 
59 template<class Type>
61 div
62 (
64 )
65 {
67  tssf.clear();
68  return Div;
69 }
70 
71 
72 template<class Type>
73 tmp
74 <
76  <
77  typename innerProduct<vector, Type>::type, fvPatchField, volMesh
78  >
79 >
80 div
81 (
83  const word& name
84 )
85 {
87  (
88  vf.mesh(), vf.mesh().schemes().div(name)
89  ).ref().fvcDiv(vf);
90 }
91 
92 
93 template<class Type>
94 tmp
95 <
97  <
98  typename innerProduct<vector, Type>::type, fvPatchField, volMesh
99  >
100 >
101 div
102 (
104  const word& name
105 )
106 {
107  typedef typename innerProduct<vector, Type>::type DivType;
109  (
110  fvc::div(tvvf(), name)
111  );
112  tvvf.clear();
113  return Div;
114 }
115 
116 template<class Type>
117 tmp
118 <
120  <
121  typename innerProduct<vector, Type>::type, fvPatchField, volMesh
122  >
123 >
124 div
125 (
127 )
128 {
129  return fvc::div(vf, "div("+vf.name()+')');
130 }
131 
132 
133 template<class Type>
134 tmp
135 <
137  <
138  typename innerProduct<vector, Type>::type, fvPatchField, volMesh
139  >
140 >
141 div
142 (
144 )
145 {
146  typedef typename innerProduct<vector, Type>::type DivType;
148  tvvf.clear();
149  return Div;
150 }
151 
152 
153 template<class Type>
155 div
156 (
157  const surfaceScalarField& flux,
159  const word& name
160 )
161 {
163  (
164  vf.mesh(),
165  flux,
166  vf.mesh().schemes().div(name)
167  ).ref().fvcDiv(flux, vf);
168 }
169 
170 
171 template<class Type>
173 div
174 (
175  const tmp<surfaceScalarField>& tflux,
177  const word& name
178 )
179 {
181  (
182  fvc::div(tflux(), vf, name)
183  );
184  tflux.clear();
185  return Div;
186 }
187 
188 
189 template<class Type>
191 div
192 (
193  const surfaceScalarField& flux,
195  const word& name
196 )
197 {
199  (
200  fvc::div(flux, tvf(), name)
201  );
202  tvf.clear();
203  return Div;
204 }
205 
206 
207 template<class Type>
209 div
210 (
211  const tmp<surfaceScalarField>& tflux,
213  const word& name
214 )
215 {
217  (
218  fvc::div(tflux(), tvf(), name)
219  );
220  tflux.clear();
221  tvf.clear();
222  return Div;
223 }
224 
225 
226 template<class Type>
228 div
229 (
230  const surfaceScalarField& flux,
232 )
233 {
234  return fvc::div
235  (
236  flux, vf, "div("+flux.name()+','+vf.name()+')'
237  );
238 }
239 
240 
241 template<class Type>
243 div
244 (
245  const tmp<surfaceScalarField>& tflux,
247 )
248 {
250  (
251  fvc::div(tflux(), vf)
252  );
253  tflux.clear();
254  return Div;
255 }
256 
257 
258 template<class Type>
260 div
261 (
262  const surfaceScalarField& flux,
264 )
265 {
267  (
268  fvc::div(flux, tvf())
269  );
270  tvf.clear();
271  return Div;
272 }
273 
274 
275 template<class Type>
277 div
278 (
279  const tmp<surfaceScalarField>& tflux,
281 )
282 {
284  (
285  fvc::div(tflux(), tvf())
286  );
287  tflux.clear();
288  tvf.clear();
289  return Div;
290 }
291 
292 
293 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
294 
295 } // End namespace fvc
296 
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
298 
299 } // End namespace Foam
300 
301 // ************************************************************************* //
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
const word & name() const
Return name.
Definition: IOobject.H:315
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
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
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
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.
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
rDeltaT ref()
Namespace for OpenFOAM.