fvmSup.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 "volFields.H"
29 #include "surfaceFields.H"
30 #include "fvMatrix.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 template<class Type>
37 (
38  const DimensionedField<Type, volMesh>& su,
39  const GeometricField<Type, fvPatchField, volMesh>& vf
40 )
41 {
42  const fvMesh& mesh = vf.mesh();
43 
44  tmp<fvMatrix<Type>> tfvm
45  (
46  new fvMatrix<Type>
47  (
48  vf,
49  dimVol*su.dimensions()
50  )
51  );
52  fvMatrix<Type>& fvm = tfvm.ref();
53 
54  fvm.source() -= mesh.V()*su.field();
55 
56  return tfvm;
57 }
58 
59 
60 template<class Type>
63 (
64  const tmp<DimensionedField<Type, volMesh>>& tsu,
65  const GeometricField<Type, fvPatchField, volMesh>& vf
66 )
67 {
68  tmp<fvMatrix<Type>> tfvm = fvm::Su(tsu(), vf);
69  tsu.clear();
70  return tfvm;
71 }
72 
73 
74 template<class Type>
77 (
78  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
79  const GeometricField<Type, fvPatchField, volMesh>& vf
80 )
81 {
82  tmp<fvMatrix<Type>> tfvm = fvm::Su(tsu(), vf);
83  tsu.clear();
84  return tfvm;
85 }
86 
87 
88 template<class Type>
91 (
92  const zero&,
93  const GeometricField<Type, fvPatchField, volMesh>& vf
94 )
95 {
96  return zeroField();
97 }
98 
99 
100 template<class Type>
103 (
104  const volScalarField::Internal& sp,
105  const GeometricField<Type, fvPatchField, volMesh>& vf
106 )
107 {
108  const fvMesh& mesh = vf.mesh();
109 
110  tmp<fvMatrix<Type>> tfvm
111  (
112  new fvMatrix<Type>
113  (
114  vf,
115  dimVol*sp.dimensions()*vf.dimensions()
116  )
117  );
118  fvMatrix<Type>& fvm = tfvm.ref();
119 
120  fvm.diag() += mesh.V()*sp.field();
121 
122  return tfvm;
123 }
124 
125 
126 template<class Type>
129 (
130  const tmp<volScalarField::Internal>& tsp,
131  const GeometricField<Type, fvPatchField, volMesh>& vf
132 )
133 {
134  tmp<fvMatrix<Type>> tfvm = fvm::Sp(tsp(), vf);
135  tsp.clear();
136  return tfvm;
137 }
138 
139 
140 template<class Type>
143 (
144  const tmp<volScalarField>& tsp,
145  const GeometricField<Type, fvPatchField, volMesh>& vf
146 )
147 {
148  tmp<fvMatrix<Type>> tfvm = fvm::Sp(tsp(), vf);
149  tsp.clear();
150  return tfvm;
151 }
152 
153 
154 template<class Type>
157 (
158  const dimensionedScalar& sp,
159  const GeometricField<Type, fvPatchField, volMesh>& vf
160 )
161 {
162  const fvMesh& mesh = vf.mesh();
163 
164  tmp<fvMatrix<Type>> tfvm
165  (
166  new fvMatrix<Type>
167  (
168  vf,
169  dimVol*sp.dimensions()*vf.dimensions()
170  )
171  );
172  fvMatrix<Type>& fvm = tfvm.ref();
173 
174  fvm.diag() += mesh.V()*sp.value();
175 
176  return tfvm;
177 }
178 
179 
180 template<class Type>
183 (
184  const zero&,
185  const GeometricField<Type, fvPatchField, volMesh>&
186 )
187 {
188  return zeroField();
189 }
190 
191 
192 template<class Type>
195 (
196  const volScalarField::Internal& susp,
197  const GeometricField<Type, fvPatchField, volMesh>& vf
198 )
199 {
200  const fvMesh& mesh = vf.mesh();
201 
202  tmp<fvMatrix<Type>> tfvm
203  (
204  new fvMatrix<Type>
205  (
206  vf,
207  dimVol*susp.dimensions()*vf.dimensions()
208  )
209  );
210  fvMatrix<Type>& fvm = tfvm.ref();
211 
212  fvm.diag() += mesh.V()*max(susp.field(), scalar(0));
213 
214  fvm.source() -= mesh.V()*min(susp.field(), scalar(0))
215  *vf.primitiveField();
216 
217  return tfvm;
218 }
219 
220 
221 template<class Type>
224 (
225  const tmp<volScalarField::Internal>& tsusp,
226  const GeometricField<Type, fvPatchField, volMesh>& vf
227 )
228 {
229  tmp<fvMatrix<Type>> tfvm = fvm::SuSp(tsusp(), vf);
230  tsusp.clear();
231  return tfvm;
232 }
233 
234 
235 template<class Type>
238 (
239  const tmp<volScalarField>& tsusp,
240  const GeometricField<Type, fvPatchField, volMesh>& vf
241 )
242 {
243  tmp<fvMatrix<Type>> tfvm = fvm::SuSp(tsusp(), vf);
244  tsusp.clear();
245  return tfvm;
246 }
247 
248 
249 template<class Type>
252 (
253  const zero&,
254  const GeometricField<Type, fvPatchField, volMesh>& vf
255 )
256 {
257  return zeroField();
258 }
259 
260 
261 // ************************************************************************* //
volFields.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::zeroField
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:47
Foam::vtk::Tools::zeroField
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Definition: foamVtkToolsTemplates.C:320
Sp
zeroField Sp
Definition: alphaSuSp.H:2
fvMatrix.H
surfaceFields.H
Foam::surfaceFields.
Foam::fvm::Su
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
Foam::fam::SuSp
tmp< faMatrix< Type > > SuSp(const areaScalarField &sp, const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famSup.C:144
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:36
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimVol
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:65