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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "sampledSurface.H"
30 
31 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
32 
33 template<class Type>
36 (
37  const interpolation<Type>& sampler,
38  const labelUList& elements,
39  const faceList& fcs,
40  const pointField& pts
41 )
42 {
43  const label len = elements.size();
44 
45  if (len != fcs.size())
46  {
48  << "size mismatch: "
49  << "sampled elements (" << len
50  << ") != faces (" << fcs.size() << ')'
51  << exit(FatalError);
52  }
53 
54  auto tvalues = tmp<Field<Type>>::New(len);
55  auto& values = tvalues.ref();
56 
57  for (label i=0; i < len; ++i)
58  {
59  const label celli = elements[i];
60  const point pt = fcs[i].centre(pts);
61 
62  values[i] = sampler.interpolate(pt, celli);
63  }
64 
65  return tvalues;
66 }
67 
68 
69 template<class Type>
72 (
74 )
75 {
76  const fvMesh& mesh = dynamic_cast<const fvMesh&>(pfld.mesh()());
77 
79  (
80  IOobject
81  (
82  "cellAvg",
83  mesh.time().timeName(),
84  pfld.db(),
85  IOobject::NO_READ,
86  IOobject::NO_WRITE,
87  false
88  ),
89  mesh,
91  );
92  auto& cellAvg = tcellAvg.ref();
93 
94  labelField nPointCells(mesh.nCells(), Zero);
95 
96  for (label pointi = 0; pointi < mesh.nPoints(); ++pointi)
97  {
98  const Type& val = pfld[pointi];
99  const labelList& pCells = mesh.pointCells(pointi);
100 
101  for (const label celli : pCells)
102  {
103  cellAvg[celli] += val;
104  ++nPointCells[celli];
105  }
106  }
107 
108  forAll(cellAvg, celli)
109  {
110  cellAvg[celli] /= nPointCells[celli];
111  }
112 
113  // Give value to calculatedFvPatchFields
114  cellAvg.correctBoundaryConditions();
115 
116  return tcellAvg;
117 }
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
122 template<class Type, class GeoMeshType>
124 (
125  const objectRegistry& obr,
126  const word& fieldName,
127  const dimensionSet& dims,
128  const Field<Type>& values,
129  word lookupName
130 ) const
131 {
132  polySurface* surfptr = this->getRegistrySurface(obr, lookupName);
133 
134  if (surfptr)
135  {
136  surfptr->storeField<Type, GeoMeshType>
137  (
138  fieldName, dims, values
139  );
140  }
141 
142  return surfptr;
143 }
144 
145 
146 template<class Type, class GeoMeshType>
148 (
149  const objectRegistry& obr,
150  const word& fieldName,
151  const dimensionSet& dims,
153  word lookupName
154 ) const
155 {
156  polySurface* surfptr = this->getRegistrySurface(obr, lookupName);
157 
158  if (surfptr)
159  {
160  surfptr->storeField<Type, GeoMeshType>
161  (
162  fieldName, dims, std::move(values)
163  );
164  }
165 
166  return surfptr;
167 }
168 
169 
170 template<class Type, class GeoMeshType>
172 (
173  const word& fieldName,
174  const dimensionSet& dims,
175  const Field<Type>& values,
176  word lookupName
177 ) const
178 {
179  surfMesh* surfptr = this->getSurfMesh(lookupName);
180 
181  if (surfptr)
182  {
183  surfptr->storeField<Type, GeoMeshType>
184  (
185  fieldName, dims, values
186  );
187  }
188 
189  return surfptr;
190 }
191 
192 
193 template<class Type, class GeoMeshType>
195 (
196  const word& fieldName,
197  const dimensionSet& dims,
199  word lookupName
200 ) const
201 {
202  surfMesh* surfptr = this->getSurfMesh(lookupName);
203 
204  if (surfptr)
205  {
206  surfptr->storeField<Type, GeoMeshType>
207  (
208  fieldName, dims, std::move(values)
209  );
210  }
211 
212  return surfptr;
213 }
214 
215 
216 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::val
label ListType::const_reference val
Definition: ListOps.H:407
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::surfMesh::storeField
void storeField(const word &fieldName, const dimensionSet &dims, const Field< Type > &values)
Copy/store named field as face or point data (template parameter).
Definition: surfMeshTemplates.C:35
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::surfMesh
A surface mesh consisting of general polygon faces.
Definition: surfMesh.H:64
Foam::polySurface::storeField
void storeField(const word &fieldName, const dimensionSet &dims, const Field< Type > &values)
Copy/store named field as face or point data (template parameter).
Definition: polySurfaceTemplates.C:96
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::sampledSurface::pointAverage
static tmp< GeometricField< Type, fvPatchField, volMesh > > pointAverage(const GeometricField< Type, pointPatchField, pointMesh > &pfld)
Create cell values by averaging the point values.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::sampledSurface::sampleOnFaces
static tmp< Field< Type > > sampleOnFaces(const interpolation< Type > &sampler, const labelUList &elements, const faceList &fcs, const pointField &pts)
General loop for sampling elements to faces.
Foam::polySurface
A surface mesh consisting of general polygon faces and capable of holding fields.
Definition: polySurface.H:67
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:62
Foam::Field< vector >
sampledSurface.H
Foam::sampledSurface::storeRegistryField
bool storeRegistryField(const objectRegistry &obr, const word &fieldName, const dimensionSet &dims, const Field< Type > &values, word lookupName="") const
Copy/store sampled field onto registered surface (if it exists)
Definition: sampledSurfaceTemplates.C:124
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:95
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Vector< scalar >
Foam::List< face >
Foam::UList< label >
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::interpolation::interpolate
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
Foam::sampledSurface::storeSurfMeshField
bool storeSurfMeshField(const word &fieldName, const dimensionSet &dims, const Field< Type > &values, word lookupName="") const
Copy/store sampled Face field onto surfMesh (if it exists)
Definition: sampledSurfaceTemplates.C:172