fvcSurfaceIntegrate.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-------------------------------------------------------------------------------
10License
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 "fvcSurfaceIntegrate.H"
29#include "fvMesh.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace fvc
40{
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44template<class Type>
46(
47 Field<Type>& ivf,
49)
50{
51 const fvMesh& mesh = ssf.mesh();
52
53 const labelUList& owner = mesh.owner();
54 const labelUList& neighbour = mesh.neighbour();
55
56 const Field<Type>& issf = ssf;
57
58 forAll(owner, facei)
59 {
60 ivf[owner[facei]] += issf[facei];
61 ivf[neighbour[facei]] -= issf[facei];
62 }
63
64 forAll(mesh.boundary(), patchi)
65 {
66 const labelUList& pFaceCells =
67 mesh.boundary()[patchi].faceCells();
68
69 const fvsPatchField<Type>& pssf = ssf.boundaryField()[patchi];
70
71 forAll(mesh.boundary()[patchi], facei)
72 {
73 ivf[pFaceCells[facei]] += pssf[facei];
74 }
75 }
76
77 ivf /= mesh.Vsc();
78}
79
80
81template<class Type>
84(
86)
87{
88 const fvMesh& mesh = ssf.mesh();
89
91 (
93 (
95 (
96 "surfaceIntegrate("+ssf.name()+')',
97 ssf.instance(),
98 mesh,
101 ),
102 mesh,
105 )
106 );
108
111
112 return tvf;
113}
114
115
116template<class Type>
119(
121)
122{
124 (
126 );
127 tssf.clear();
128 return tvf;
129}
130
131
132template<class Type>
135(
137)
138{
139 const fvMesh& mesh = ssf.mesh();
140
142 (
144 (
146 (
147 "surfaceSum("+ssf.name()+')',
148 ssf.instance(),
149 mesh,
152 ),
153 mesh,
156 )
157 );
159
160 const labelUList& owner = mesh.owner();
161 const labelUList& neighbour = mesh.neighbour();
162
163 forAll(owner, facei)
164 {
165 vf[owner[facei]] += ssf[facei];
166 vf[neighbour[facei]] += ssf[facei];
167 }
168
169 forAll(mesh.boundary(), patchi)
170 {
171 const labelUList& pFaceCells =
172 mesh.boundary()[patchi].faceCells();
173
174 const fvsPatchField<Type>& pssf = ssf.boundaryField()[patchi];
175
176 forAll(mesh.boundary()[patchi], facei)
177 {
178 vf[pFaceCells[facei]] += pssf[facei];
179 }
180 }
181
183
184 return tvf;
185}
186
187
188template<class Type>
190(
192)
193{
195 tssf.clear();
196 return tvf;
197}
198
199
200// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201
202} // End namespace fvc
203
204// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205
206} // End namespace Foam
207
208// ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
Generic dimensioned Type class.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
dynamicFvMesh & mesh
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:64
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333