sampledSurfaceTemplates.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-2015 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 "sampledSurface.H"
27 
28 template<class Type>
30 {
31  if (faces().empty() || field.empty())
32  {
33  return false;
34  }
35 
36  if (field.size() != faces().size())
37  {
39  << "size mismatch: "
40  << "field (" << field.size()
41  << ") != surface (" << faces().size() << ")"
42  << exit(FatalError);
43  }
44 
45  return true;
46 }
47 
48 
49 template<class Type>
51 {
52  Type value = pTraits<Type>::zero;
53 
54  if (checkFieldSize(field))
55  {
56  value = sum(field*magSf());
57  }
58 
59  reduce(value, sumOp<Type>());
60  return value;
61 }
62 
63 
64 template<class Type>
66 {
67  Type value = integrate(field());
68  field.clear();
69  return value;
70 }
71 
72 
73 template<class Type>
75 {
76  Type value = pTraits<Type>::zero;
77 
78  if (checkFieldSize(field))
79  {
80  value = sum(field*magSf());
81  }
82 
83  reduce(value, sumOp<Type>());
84 
85  // avoid divide-by-zero
86  if (area())
87  {
88  return value/area();
89  }
90  else
91  {
92  return pTraits<Type>::zero;
93  }
94 }
95 
96 
97 template<class Type>
98 Type Foam::sampledSurface::average(const tmp<Field<Type> >& field) const
99 {
100  Type value = average(field());
101  field.clear();
102  return value;
103 }
104 
105 
106 template<class ReturnType, class Type>
108 (
109  Field<ReturnType>& res,
110  const Field<Type>& field
111 ) const
112 {
113  if (checkFieldSize(field))
114  {
115  const vectorField& norm = Sf();
116 
117  forAll(norm, faceI)
118  {
119  res[faceI] = field[faceI] & (norm[faceI]/mag(norm[faceI]));
120  }
121  }
122  else
123  {
124  res.clear();
125  }
126 }
127 
128 
129 template<class ReturnType, class Type>
131 (
132  Field<ReturnType>& res,
133  const tmp<Field<Type> >& field
134 ) const
135 {
136  project(res, field());
137  field.clear();
138 }
139 
140 
141 template<class ReturnType, class Type>
144 (
145  const tmp<Field<Type> >& field
146 ) const
147 {
148  tmp<Field<ReturnType> > tRes(new Field<ReturnType>(faces().size()));
149  project(tRes(), field);
150  return tRes;
151 }
152 
153 
154 template<class Type>
157 (
159 ) const
160 {
161  const fvMesh& mesh = dynamic_cast<const fvMesh&>(pfld.mesh()());
162 
164  (
166  (
167  IOobject
168  (
169  "cellAvg",
170  mesh.time().timeName(),
171  pfld.db(),
174  false
175  ),
176  mesh,
178  )
179  );
180  GeometricField<Type, fvPatchField, volMesh>& cellAvg = tcellAvg();
181 
182  labelField nPointCells(mesh.nCells(), 0);
183  {
184  for (label pointI = 0; pointI < mesh.nPoints(); pointI++)
185  {
186  const labelList& pCells = mesh.pointCells(pointI);
187 
188  forAll(pCells, i)
189  {
190  label cellI = pCells[i];
191 
192  cellAvg[cellI] += pfld[pointI];
193  nPointCells[cellI]++;
194  }
195  }
196  }
197  forAll(cellAvg, cellI)
198  {
199  cellAvg[cellI] /= nPointCells[cellI];
200  }
201  // Give value to calculatedFvPatchFields
202  cellAvg.correctBoundaryConditions();
203 
204  return tcellAvg;
205 }
206 
207 
208 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::sampledSurface::faces
virtual const faceList & faces() const =0
Faces of surface.
Foam::sampledSurface::pointAverage
tmp< GeometricField< Type, fvPatchField, volMesh > > pointAverage(const GeometricField< Type, pointPatchField, pointMesh > &pfld) const
Interpolate from points to cell centre.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::sampledSurface::average
Type average(const Field< Type > &) const
Area-averaged value of a field across the surface.
Definition: sampledSurfaceTemplates.C:74
Foam::sampledSurface::project
void project(Field< ReturnType > &, const Field< Type > &) const
Project field onto surface.
Definition: sampledSurfaceTemplates.C:108
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::sampledSurface::integrate
Type integrate(const Field< Type > &) const
Integration of a field across the surface.
Definition: sampledSurfaceTemplates.C:50
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field< Type >
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
sampledSurface.H
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned< Type >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::primitiveMesh::pointCells
const labelListList & pointCells() const
Definition: primitiveMeshPointCells.C:108
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::sampledSurface::checkFieldSize
bool checkFieldSize(const Field< Type > &) const
Check field size matches surface size.
Definition: sampledSurfaceTemplates.C:29
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:335