fvcAverage.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-2019 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 "fvcAverage.H"
27 #include "fvcSurfaceIntegrate.H"
28 #include "fvMesh.H"
29 #include "linear.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace fvc
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 template<class Type>
44 tmp<GeometricField<Type, fvPatchField, volMesh>>
45 average
46 (
48 )
49 {
50  const fvMesh& mesh = ssf.mesh();
51 
53  (
55  (
56  "average("+ssf.name()+')',
57  mesh,
58  dimensioned<Type>("0", ssf.dimensions(), Zero)
59  )
60  );
61 
62  if (!mesh.nGeometricD())
63  {
64  return taverage;
65  }
66 
68 
69  av.primitiveFieldRef() =
70  (
71  surfaceSum(mesh.magSf()*ssf)().primitiveField()
72  /surfaceSum(mesh.magSf())().primitiveField()
73  );
74 
76  Boundary& bav = av.boundaryFieldRef();
77 
78  forAll(bav, patchi)
79  {
80  bav[patchi] = ssf.boundaryField()[patchi];
81  }
82 
84 
85  return taverage;
86 }
87 
88 
89 template<class Type>
91 average
92 (
94 )
95 {
97  (
98  fvc::average(tssf())
99  );
100  tssf.clear();
101  return taverage;
102 }
103 
104 
105 template<class Type>
107 average
108 (
110 )
111 {
112  return fvc::average(linearInterpolate(vtf));
113 }
114 
115 
116 template<class Type>
118 average
119 (
121 )
122 {
124  (
125  fvc::average(tvtf())
126  );
127  tvtf.clear();
128  return taverage;
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 } // End namespace fvc
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 } // End namespace Foam
139 
140 // ************************************************************************* //
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const word & name() const
Return name.
Definition: IOobject.H:303
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:108
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Generic GeometricField class.
Generic dimensioned Type class.
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:825
Area-weighted average a surfaceField creating a volField.
const dimensionSet & dimensions() const
Return dimensions.
dynamicFvMesh & mesh
static const zero Zero
Definition: zero.H:97
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
const Mesh & mesh() const
Return mesh.
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void correctBoundaryConditions()
Correct boundary field.
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.