fvcAverage.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 | Copyright (C) 2011-2016 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  IOobject
57  (
58  "average("+ssf.name()+')',
59  ssf.instance(),
60  mesh,
63  ),
64  mesh,
65  ssf.dimensions()
66  )
67  );
69 
70  av.primitiveFieldRef() =
71  (
72  surfaceSum(mesh.magSf()*ssf)().primitiveField()
73  /surfaceSum(mesh.magSf())().primitiveField()
74  );
75 
77  Boundary& bav = av.boundaryFieldRef();
78 
79  forAll(bav, patchi)
80  {
81  bav[patchi] = ssf.boundaryField()[patchi];
82  }
83 
85 
86  return taverage;
87 }
88 
89 
90 template<class Type>
92 average
93 (
95 )
96 {
98  (
99  fvc::average(tssf())
100  );
101  tssf.clear();
102  return taverage;
103 }
104 
105 
106 template<class Type>
108 average
109 (
111 )
112 {
113  return fvc::average(linearInterpolate(vtf));
114 }
115 
116 
117 template<class Type>
119 average
120 (
122 )
123 {
125  (
126  fvc::average(tvtf())
127  );
128  tvtf.clear();
129  return taverage;
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 
135 } // End namespace fvc
136 
137 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138 
139 } // End namespace Foam
140 
141 // ************************************************************************* //
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:428
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:107
Generic GeometricField class.
Area-weighted average a surfaceField creating a volField.
dynamicFvMesh & mesh
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
const dimensionSet & dimensions() const
Return dimensions.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const Mesh & mesh() const
Return mesh.
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:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const fileName & instance() const
Definition: IOobject.H:337
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.