All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
fvmSup.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 "volFields.H"
27 #include "surfaceFields.H"
28 #include "fvMatrix.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 template<class Type>
35 (
36  const DimensionedField<Type, volMesh>& su,
37  const GeometricField<Type, fvPatchField, volMesh>& vf
38 )
39 {
40  const fvMesh& mesh = vf.mesh();
41 
42  tmp<fvMatrix<Type>> tfvm
43  (
44  new fvMatrix<Type>
45  (
46  vf,
47  dimVol*su.dimensions()
48  )
49  );
50  fvMatrix<Type>& fvm = tfvm.ref();
51 
52  fvm.source() -= mesh.V()*su.field();
53 
54  return tfvm;
55 }
56 
57 
58 template<class Type>
61 (
62  const tmp<DimensionedField<Type, volMesh>>& tsu,
63  const GeometricField<Type, fvPatchField, volMesh>& vf
64 )
65 {
66  tmp<fvMatrix<Type>> tfvm = fvm::Su(tsu(), vf);
67  tsu.clear();
68  return tfvm;
69 }
70 
71 
72 template<class Type>
75 (
76  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
77  const GeometricField<Type, fvPatchField, volMesh>& vf
78 )
79 {
80  tmp<fvMatrix<Type>> tfvm = fvm::Su(tsu(), vf);
81  tsu.clear();
82  return tfvm;
83 }
84 
85 
86 template<class Type>
89 (
90  const zero&,
91  const GeometricField<Type, fvPatchField, volMesh>& vf
92 )
93 {
94  return zeroField();
95 }
96 
97 
98 template<class Type>
101 (
102  const volScalarField::Internal& sp,
103  const GeometricField<Type, fvPatchField, volMesh>& vf
104 )
105 {
106  const fvMesh& mesh = vf.mesh();
107 
108  tmp<fvMatrix<Type>> tfvm
109  (
110  new fvMatrix<Type>
111  (
112  vf,
113  dimVol*sp.dimensions()*vf.dimensions()
114  )
115  );
116  fvMatrix<Type>& fvm = tfvm.ref();
117 
118  fvm.diag() += mesh.V()*sp.field();
119 
120  return tfvm;
121 }
122 
123 
124 template<class Type>
127 (
128  const tmp<volScalarField::Internal>& tsp,
129  const GeometricField<Type, fvPatchField, volMesh>& vf
130 )
131 {
132  tmp<fvMatrix<Type>> tfvm = fvm::Sp(tsp(), vf);
133  tsp.clear();
134  return tfvm;
135 }
136 
137 
138 template<class Type>
141 (
142  const tmp<volScalarField>& tsp,
143  const GeometricField<Type, fvPatchField, volMesh>& vf
144 )
145 {
146  tmp<fvMatrix<Type>> tfvm = fvm::Sp(tsp(), vf);
147  tsp.clear();
148  return tfvm;
149 }
150 
151 
152 template<class Type>
155 (
156  const dimensionedScalar& sp,
157  const GeometricField<Type, fvPatchField, volMesh>& vf
158 )
159 {
160  const fvMesh& mesh = vf.mesh();
161 
162  tmp<fvMatrix<Type>> tfvm
163  (
164  new fvMatrix<Type>
165  (
166  vf,
167  dimVol*sp.dimensions()*vf.dimensions()
168  )
169  );
170  fvMatrix<Type>& fvm = tfvm.ref();
171 
172  fvm.diag() += mesh.V()*sp.value();
173 
174  return tfvm;
175 }
176 
177 
178 template<class Type>
181 (
182  const zero&,
183  const GeometricField<Type, fvPatchField, volMesh>&
184 )
185 {
186  return zeroField();
187 }
188 
189 
190 template<class Type>
193 (
194  const volScalarField::Internal& susp,
195  const GeometricField<Type, fvPatchField, volMesh>& vf
196 )
197 {
198  const fvMesh& mesh = vf.mesh();
199 
200  tmp<fvMatrix<Type>> tfvm
201  (
202  new fvMatrix<Type>
203  (
204  vf,
205  dimVol*susp.dimensions()*vf.dimensions()
206  )
207  );
208  fvMatrix<Type>& fvm = tfvm.ref();
209 
210  fvm.diag() += mesh.V()*max(susp.field(), scalar(0));
211 
212  fvm.source() -= mesh.V()*min(susp.field(), scalar(0))
213  *vf.primitiveField();
214 
215  return tfvm;
216 }
217 
218 
219 template<class Type>
222 (
223  const tmp<volScalarField::Internal>& tsusp,
224  const GeometricField<Type, fvPatchField, volMesh>& vf
225 )
226 {
227  tmp<fvMatrix<Type>> tfvm = fvm::SuSp(tsusp(), vf);
228  tsusp.clear();
229  return tfvm;
230 }
231 
232 
233 template<class Type>
236 (
237  const tmp<volScalarField>& tsusp,
238  const GeometricField<Type, fvPatchField, volMesh>& vf
239 )
240 {
241  tmp<fvMatrix<Type>> tfvm = fvm::SuSp(tsusp(), vf);
242  tsusp.clear();
243  return tfvm;
244 }
245 
246 
247 template<class Type>
250 (
251  const zero&,
252  const GeometricField<Type, fvPatchField, volMesh>& vf
253 )
254 {
255  return zeroField();
256 }
257 
258 
259 // ************************************************************************* //
Foam::surfaceFields.
tmp< GeometricField< Type, fvPatchField, volMesh > > SuSp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:102
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
zeroField Su
Definition: alphaSuSp.H:1
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:50
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
zeroField Sp
Definition: alphaSuSp.H:2