sampledSurfaceTemplates.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-2022 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Type>
31 bool Foam::sampledSurface::checkFieldSize(const Field<Type>& field) const
32 {
33  if (faces().empty() || field.empty())
34  {
35  return false;
36  }
37 
38  if (field.size() != faces().size())
39  {
41  << "size mismatch: "
42  << "field (" << field.size()
43  << ") != surface (" << faces().size() << ")"
44  << exit(FatalError);
45  }
46 
47  return true;
48 }
49 
50 
51 template<class Type>
53 {
54  Type value = Zero;
55 
56  if (checkFieldSize(field))
57  {
58  value = sum(field*magSf());
59  }
60 
61  reduce(value, sumOp<Type>());
62  return value;
63 }
64 
65 
66 template<class Type>
68 {
69  Type value = integrate(field());
70  field.clear();
71  return value;
72 }
73 
74 
75 template<class Type>
77 {
78  Type value = Zero;
79 
80  if (checkFieldSize(field))
81  {
82  value = sum(field*magSf());
83  }
84 
85  reduce(value, sumOp<Type>());
86 
87  // avoid divide-by-zero
88  if (area())
89  {
90  return value/area();
91  }
92  else
93  {
94  return Zero;
95  }
96 }
97 
98 
99 template<class Type>
101 {
102  Type value = average(field());
103  field.clear();
104  return value;
105 }
106 
107 
108 template<class ReturnType, class Type>
109 void Foam::sampledSurface::project
110 (
111  Field<ReturnType>& res,
112  const Field<Type>& field
113 ) const
114 {
115  if (checkFieldSize(field))
116  {
117  const vectorField& norm = Sf();
118 
119  forAll(norm, facei)
120  {
121  res[facei] = field[facei] & (norm[facei]/mag(norm[facei]));
122  }
123  }
124  else
125  {
126  res.clear();
127  }
128 }
129 
130 
131 template<class ReturnType, class Type>
132 void Foam::sampledSurface::project
133 (
134  Field<ReturnType>& res,
135  const tmp<Field<Type>>& field
136 ) const
137 {
138  project(res, field());
139  field.clear();
140 }
141 
142 
143 template<class ReturnType, class Type>
145 Foam::sampledSurface::project
146 (
147  const tmp<Field<Type>>& field
148 ) const
149 {
150  tmp<Field<ReturnType>> tRes(new Field<ReturnType>(faces().size()));
151  project(tRes(), field);
152  return tRes;
153 }
154 
155 
156 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Pre-declare SubField and related Field type.
Definition: Field.H:82
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
virtual const faceList & faces() const =0
Faces of surface.
Type average(const Field< Type > &) const
Area-averaged value of a field across the surface.
Type integrate(const Field< Type > &) const
Integration of a field across the surface.
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError