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 Copyright (C) 2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "volFields.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33template<class Type>
36(
37 const Foam::zero,
38 const GeometricField<Type, fvPatchField, volMesh>& fld
39)
40{
41 return zeroField();
42}
43
44
45template<class Type>
48(
49 const dimensioned<Type>& su,
50 const GeometricField<Type, fvPatchField, volMesh>& fld
51)
52{
53 auto tmat = tmp<fvMatrix<Type>>::New
54 (
55 fld,
56 dimVol*su.dimensions()
57 );
58 auto& mat = tmat.ref();
59 const auto& domain = fld.mesh().V();
60
61 if (magSqr(su.value()) > VSMALL)
62 {
63 mat.source() -= domain*su.value();
64 }
65
66 return tmat;
67}
68
69
70template<class Type>
73(
74 const DimensionedField<Type, volMesh>& su,
75 const GeometricField<Type, fvPatchField, volMesh>& fld
76)
77{
78 auto tmat = tmp<fvMatrix<Type>>::New
79 (
80 fld,
81 dimVol*su.dimensions()
82 );
83 auto& mat = tmat.ref();
84 const auto& domain = fld.mesh().V();
85
86 mat.source() -= domain*su.field();
87
88 return tmat;
89}
90
91
92template<class Type>
95(
96 const tmp<DimensionedField<Type, volMesh>>& tsu,
97 const GeometricField<Type, fvPatchField, volMesh>& fld
98)
99{
100 tmp<fvMatrix<Type>> tmat = fvm::Su(tsu(), fld);
101 tsu.clear();
102 return tmat;
103}
104
105
106template<class Type>
109(
110 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
111 const GeometricField<Type, fvPatchField, volMesh>& fld
112)
113{
114 tmp<fvMatrix<Type>> tmat = fvm::Su(tsu(), fld);
115 tsu.clear();
116 return tmat;
117}
118
119
120// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121
122template<class Type>
125(
126 const Foam::zero,
127 const GeometricField<Type, fvPatchField, volMesh>&
128)
129{
130 return zeroField();
131}
132
133
134template<class Type>
137(
138 const dimensionedScalar& sp,
139 const GeometricField<Type, fvPatchField, volMesh>& fld
140)
141{
142 auto tmat = tmp<fvMatrix<Type>>::New
143 (
144 fld,
145 dimVol*sp.dimensions()*fld.dimensions()
146 );
147 auto& mat = tmat.ref();
148 const auto& domain = fld.mesh().V();
149
150 if (mag(sp.value()) > ROOTVSMALL)
151 {
152 mat.diag() += domain*sp.value();
153 }
154
155 return tmat;
156}
157
158
159template<class Type>
162(
163 const DimensionedField<scalar, volMesh>& sp,
164 const GeometricField<Type, fvPatchField, volMesh>& fld
165)
166{
167 auto tmat = tmp<fvMatrix<Type>>::New
168 (
169 fld,
170 dimVol*sp.dimensions()*fld.dimensions()
171 );
172 auto& mat = tmat.ref();
173 const auto& domain = fld.mesh().V();
174
175 mat.diag() += domain*sp.field();
176
177 return tmat;
178}
179
180
181template<class Type>
184(
185 const tmp<DimensionedField<scalar, volMesh>>& tsp,
186 const GeometricField<Type, fvPatchField, volMesh>& fld
187)
188{
189 tmp<fvMatrix<Type>> tmat = fvm::Sp(tsp(), fld);
190 tsp.clear();
191 return tmat;
192}
193
194
195template<class Type>
198(
199 const tmp<volScalarField>& tsp,
200 const GeometricField<Type, fvPatchField, volMesh>& fld
201)
202{
203 tmp<fvMatrix<Type>> tmat = fvm::Sp(tsp(), fld);
204 tsp.clear();
205 return tmat;
206}
207
208
209// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210
211template<class Type>
214(
215 const Foam::zero,
216 const GeometricField<Type, fvPatchField, volMesh>& fld
217)
218{
219 return zeroField();
220}
221
222
223template<class Type>
226(
227 const dimensionedScalar& susp,
228 const GeometricField<Type, fvPatchField, volMesh>& fld
229)
230{
231 auto tmat = tmp<fvMatrix<Type>>::New
232 (
233 fld,
234 dimVol*susp.dimensions()*fld.dimensions()
235 );
236 auto& mat = tmat.ref();
237 const auto& domain = fld.mesh().V();
238
239 if (susp.value() > ROOTVSMALL)
240 {
241 mat.diag() += domain*susp.value();
242 }
243 else if (susp.value() < -ROOTVSMALL)
244 {
245 mat.source() -= domain*susp.value()*fld.primitiveField();
246 }
247
248 return tmat;
249}
250
251
252template<class Type>
255(
256 const DimensionedField<scalar, volMesh>& susp,
257 const GeometricField<Type, fvPatchField, volMesh>& fld
258)
259{
260 auto tmat = tmp<fvMatrix<Type>>::New
261 (
262 fld,
263 dimVol*susp.dimensions()*fld.dimensions()
264 );
265 auto& mat = tmat.ref();
266 const auto& domain = fld.mesh().V();
267
268 mat.diag() += domain*max(susp.field(), scalar(0));
269
270 mat.source() -= domain*min(susp.field(), scalar(0))*fld.primitiveField();
271
272 return tmat;
273}
274
275
276template<class Type>
279(
280 const tmp<DimensionedField<scalar, volMesh>>& tsusp,
281 const GeometricField<Type, fvPatchField, volMesh>& fld
282)
283{
284 tmp<fvMatrix<Type>> tmat = fvm::SuSp(tsusp(), fld);
285 tsusp.clear();
286 return tmat;
287}
288
289
290template<class Type>
293(
294 const tmp<volScalarField>& tsusp,
295 const GeometricField<Type, fvPatchField, volMesh>& fld
296)
297{
298 tmp<fvMatrix<Type>> tmat = fvm::SuSp(tsusp(), fld);
299 tsusp.clear();
300 return tmat;
301}
302
303
304// ************************************************************************* //
Y[inertIndex] max(0.0)
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A class for managing temporary objects.
Definition: tmp.H:65
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:56
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
zeroField Su(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)