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<VolField<Type>>
47 (
48  const SurfaceField<Type>& ssf
49 )
50 {
51  return VolField<Type>::New
52  (
53  "div("+ssf.name()+')',
55  );
56 }
57 
58 
59 template<class Type>
62 (
63  const tmp<SurfaceField<Type>>& tssf
64 )
65 {
66  tmp<VolField<Type>> Div(fvc::div(tssf()));
67  tssf.clear();
68  return Div;
69 }
70 
71 
72 template<class Type>
75 (
76  const VolField<Type>& vf,
77  const word& name
78 )
79 {
81  (
82  vf.mesh(), vf.mesh().schemes().div(name)
83  ).ref().fvcDiv(vf);
84 }
85 
86 
87 template<class Type>
90 (
91  const tmp<VolField<Type>>& tvvf,
92  const word& name
93 )
94 {
95  typedef typename innerProduct<vector, Type>::type DivType;
97  (
98  fvc::div(tvvf(), name)
99  );
100  tvvf.clear();
101  return Div;
102 }
103 
104 template<class Type>
107 (
108  const VolField<Type>& vf
109 )
110 {
111  return fvc::div(vf, "div("+vf.name()+')');
112 }
113 
114 
115 template<class Type>
118 (
119  const tmp<VolField<Type>>& tvvf
120 )
121 {
122  typedef typename innerProduct<vector, Type>::type DivType;
123  tmp<VolField<DivType>> Div(fvc::div(tvvf()));
124  tvvf.clear();
125  return Div;
126 }
127 
128 
129 template<class Type>
132 (
133  const surfaceScalarField& flux,
134  const VolField<Type>& vf,
135  const word& name
136 )
137 {
139  (
140  vf.mesh(),
141  flux,
142  vf.mesh().schemes().div(name)
143  ).ref().fvcDiv(flux, vf);
144 }
145 
146 
147 template<class Type>
150 (
151  const tmp<surfaceScalarField>& tflux,
152  const VolField<Type>& vf,
153  const word& name
154 )
155 {
156  tmp<VolField<Type>> Div
157  (
158  fvc::div(tflux(), vf, name)
159  );
160  tflux.clear();
161  return Div;
162 }
163 
164 
165 template<class Type>
168 (
169  const surfaceScalarField& flux,
170  const tmp<VolField<Type>>& tvf,
171  const word& name
172 )
173 {
174  tmp<VolField<Type>> Div
175  (
176  fvc::div(flux, tvf(), name)
177  );
178  tvf.clear();
179  return Div;
180 }
181 
182 
183 template<class Type>
186 (
187  const tmp<surfaceScalarField>& tflux,
188  const tmp<VolField<Type>>& tvf,
189  const word& name
190 )
191 {
192  tmp<VolField<Type>> Div
193  (
194  fvc::div(tflux(), tvf(), name)
195  );
196  tflux.clear();
197  tvf.clear();
198  return Div;
199 }
200 
201 
202 template<class Type>
205 (
206  const surfaceScalarField& flux,
207  const VolField<Type>& vf
208 )
209 {
210  return fvc::div
211  (
212  flux, vf, "div("+flux.name()+','+vf.name()+')'
213  );
214 }
215 
216 
217 template<class Type>
220 (
221  const tmp<surfaceScalarField>& tflux,
222  const VolField<Type>& vf
223 )
224 {
225  tmp<VolField<Type>> Div
226  (
227  fvc::div(tflux(), vf)
228  );
229  tflux.clear();
230  return Div;
231 }
232 
233 
234 template<class Type>
237 (
238  const surfaceScalarField& flux,
239  const tmp<VolField<Type>>& tvf
240 )
241 {
242  tmp<VolField<Type>> Div
243  (
244  fvc::div(flux, tvf())
245  );
246  tvf.clear();
247  return Div;
248 }
249 
250 
251 template<class Type>
254 (
255  const tmp<surfaceScalarField>& tflux,
256  const tmp<VolField<Type>>& tvf
257 )
258 {
259  tmp<VolField<Type>> Div
260  (
261  fvc::div(tflux(), tvf())
262  );
263  tflux.clear();
264  tvf.clear();
265  return Div;
266 }
267 
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 } // End namespace fvc
272 
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274 
275 } // End namespace Foam
276 
277 // ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
const word & name() const
Return name.
Definition: IOobject.H:310
static tmp< convectionScheme< Type > > New(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
Return a pointer to a new convectionScheme created on freestore.
static tmp< divScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new divScheme created on freestore.
Definition: divScheme.C:47
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
A class for handling words, derived from string.
Definition: word.H:62
Calculate the divergence of the given field.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
void surfaceIntegrate(Field< Type > &ivf, const SurfaceField< Type > &ssf)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
Namespace for OpenFOAM.
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39