sampledMeshedSurface.H
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) 2016-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::sampledMeshedSurface
29
30Description
31 A sampledSurface from a meshed surface.
32 It samples on the points/faces of the meshed surface.
33
34 - it either samples cells or (non-coupled) boundary faces
35
36 - 6 different modes:
37 - source=cells, interpolate=false:
38 finds per surface face centre the \em nearest cell centre
39 and uses its value
40 - source=cells, interpolate=true:
41 finds per surface face centre the \em nearest cell centre.
42 Per surface point checks if this nearest cell is the one containing
43 point; otherwise projects the point onto the nearest point on
44 the boundary of the cell (to make sure interpolateCellPoint
45 gets a valid location)
46
47 - source=insideCells, interpolate=false:
48 finds per surface face centre the cell containing it
49 and uses its value.
50 Trims surface faces outside of the mesh.
51 - source=insideCells, interpolate=true:
52 Per surface point interpolate cell containing it.
53
54 - source=boundaryFaces, interpolate=false:
55 finds per surface face centre the \em nearest point on the boundary
56 (uncoupled faces only) and uses the value (or 0 if the nearest
57 is on an empty boundary)
58 - source=boundaryFaces, interpolate=true:
59 finds per surface face centre the \em nearest point on the boundary
60 (uncoupled faces only).
61 Per surface point projects the point onto this boundary face
62 (to make sure interpolateCellPoint gets a valid location)
63
64 - since it finds the nearest per surface face, each surface face
65 is guaranteed to be on one processor only.
66 So after stitching the original surface should be complete.
67
68 This is often embedded as part of a sampled surfaces function object.
69
70Usage
71 Example of function object partial specification:
72 \verbatim
73 surfaces
74 {
75 surface1
76 {
77 type meshedSurface;
78 surface something.obj;
79 source cells;
80
81 //- Max sampling distance
82 maxDistance 0.005;
83
84 //- Fallback for missed sampling in 'cells' mode
85 defaultValue
86 {
87 "p.*" 1e5;
88 T 273.15;
89 }
90 }
91 }
92 \endverbatim
93
94 Where the sub-entries comprise:
95 \table
96 Property | Description | Required | Default
97 type | meshedSurface | yes |
98 surface | surface name in triSurface/ | yes |
99 patches | Limit to named surface regions (wordRes) | no |
100 source | cells/insideCells/boundaryFaces | yes |
101 keepIds | pass through id numbering | no | true
102 file | Alternative file name | no |
103 fileType | The surface format | no | (extension)
104 scale | Surface scaling factor | no | 0
105 maxDistance | Max search distance | no | GREAT
106 defaultValue | Value beyond max distance (dictionary) | no | empty
107 \endtable
108
109
110SourceFiles
111 sampledMeshedSurface.C
112 sampledMeshedSurfaceTemplates.C
113
114\*---------------------------------------------------------------------------*/
115
116#ifndef sampledMeshedSurface_H
117#define sampledMeshedSurface_H
118
119#include "sampledSurface.H"
120#include "MeshedSurfaces.H"
121
122// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123
124namespace Foam
125{
126
127// Forward Declarations
128class meshSearch;
129
130/*---------------------------------------------------------------------------*\
131 Class sampledMeshedSurface Declaration
132\*---------------------------------------------------------------------------*/
133
134class sampledMeshedSurface
135:
136 public sampledSurface,
137 public meshedSurface
138{
139public:
140
141 //- Types of sampling regions
142 enum class samplingSource : unsigned char
143 {
144 cells,
147 };
148
149private:
150
151 //- The concrete storage type for faces and points
152 typedef meshedSurface Mesh;
153
154
155 // Private Data
156
157 static const Enum<samplingSource> samplingSourceNames_;
158
159 //- The name of the input surface
160 word surfaceName_;
161
162 //- Surface to sample on
163 meshedSurface surface_;
164
165 //- Whether to sample internal cell values or boundary values
166 const samplingSource sampleSource_;
167
168 //- Track if the surface needs an update
169 mutable bool needsUpdate_;
170
171 //- Retain element ids/order of original surface
172 bool keepIds_;
173
174 //- For compatibility with the meshSurf interface
175 labelList zoneIds_;
176
177 //- From local surface triangle to mesh cell/face.
178 labelList sampleElements_;
179
180 //- Local points to sample per point
181 pointField samplePoints_;
182
183 //- Max search distance squared
184 scalar maxDistanceSqr_;
185
186 //- Out-of-range fallback values (when beyond the max search distance)
187 dictionary defaultValues_;
189
190 // Private Member Functions
191
192 //- Set zoneIds list based on the surfZone information
193 void setZoneMap();
194
195 //- Sample volume field onto surface faces
196 template<class Type>
197 tmp<Field<Type>> sampleOnFaces
198 (
199 const interpolation<Type>& sampler
200 ) const;
201
202 //- Interpolate volume field onto surface points
203 template<class Type>
204 tmp<Field<Type>> sampleOnPoints
205 (
206 const interpolation<Type>& interpolator
207 ) const;
208
209 bool update(const meshSearch& meshSearcher);
210
211public:
212
213 //- Declare type-name, virtual type (with debug switch)
214 TypeName("sampledMeshedSurface");
215
216
217 // Constructors
218
219 //- Construct from components
221 (
222 const word& name,
223 const polyMesh& mesh,
224 const word& surfaceName,
225 const samplingSource sampleSource
226 );
227
228 //- Construct from dictionary
230 (
231 const word& name,
232 const polyMesh& mesh,
233 const dictionary& dict
234 );
235
236
237 //- Destructor
238 virtual ~sampledMeshedSurface() = default;
239
240
241 // Member Functions
242
243 //- Does the surface need an update?
244 virtual bool needsUpdate() const;
245
246 //- Mark the surface as needing an update.
247 // May also free up unneeded data.
248 // Return false if surface was already marked as expired.
249 virtual bool expire();
250
251 //- Update the surface as required.
252 // Do nothing (and return false) if no update was needed
253 virtual bool update();
254
255 //- Update the surface using a bound box to limit the searching.
256 // For direct use, i.e. not through sample.
257 // Do nothing (and return false) if no update was needed
258 bool update(const treeBoundBox& bb);
259
260 //- Points of surface
261 virtual const pointField& points() const
262 {
263 return Mesh::points();
264 }
265
266 //- Faces of surface
267 virtual const faceList& faces() const
269 return Mesh::surfFaces();
270 }
271
272 //- Per-face zone/region information
273 virtual const labelList& zoneIds() const
274 {
275 return zoneIds_;
276 }
277
278 //- Per-face identifier (eg, element Id)
279 virtual const labelList& faceIds() const
280 {
281 return Mesh::faceIds();
282 }
283
284 //- Face area vectors
285 virtual const vectorField& Sf() const
286 {
287 return Mesh::Sf();
288 }
289
290 //- Face area magnitudes
291 virtual const scalarField& magSf() const
293 return Mesh::magSf();
294 }
295
296 //- Face centres
297 virtual const vectorField& Cf() const
298 {
299 return Mesh::Cf();
300 }
301
302 //- If element ids/order of the original surface are kept
303 virtual bool hasFaceIds() const
304 {
305 return keepIds_ && !Mesh::faceIds().empty();
306 }
307
308 //- Sampling boundary values instead of cell values
309 bool onBoundary() const
310 {
311 return sampleSource_ == samplingSource::boundaryFaces;
312 }
313
314
315 // Sample
316
317 //- Sample volume field onto surface faces
319 (
320 const interpolation<scalar>& sampler
321 ) const;
322
323 //- Sample volume field onto surface faces
325 (
326 const interpolation<vector>& sampler
327 ) const;
328
329 //- Sample volume field onto surface faces
331 (
332 const interpolation<sphericalTensor>& sampler
333 ) const;
334
335 //- Sample volume field onto surface faces
337 (
338 const interpolation<symmTensor>& sampler
339 ) const;
340
341 //- Sample volume field onto surface faces
343 (
344 const interpolation<tensor>& sampler
345 ) const;
346
347
348 // Interpolate
349
350 //- Interpolate volume field onto surface points
352 (
353 const interpolation<scalar>& interpolator
354 ) const;
355
356 //- Interpolate volume field onto surface points
358 (
359 const interpolation<vector>& interpolator
360 ) const;
361
362 //- Interpolate volume field onto surface points
364 (
365 const interpolation<sphericalTensor>& interpolator
366 ) const;
367
368 //- Interpolate volume field onto surface points
370 (
371 const interpolation<symmTensor>& interpolator
372 ) const;
373
374 //- Interpolate volume field onto surface points
376 (
377 const interpolation<tensor>& interpolator
378 ) const;
379
380
381 // Output
382
383 //- Print information
384 virtual void print(Ostream& os, int level=0) const;
385};
386
387
388// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
389
390} // End namespace Foam
391
392// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
393
394#ifdef NoRepository
396#endif
397
398// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399
400#endif
401
402// ************************************************************************* //
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
Minimal example by using system/controlDict.functions:
const vectorField & Sf() const
Face area vectors (normals)
const scalarField & magSf() const
Face area magnitudes.
const vectorField & Cf() const
Face centres.
const labelList & faceIds() const
Return const access to faces ids.
const List< face > & surfFaces() const
Return const access to the faces.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A sampledSurface from a meshed surface. It samples on the points/faces of the meshed surface.
virtual const pointField & points() const
Points of surface.
bool onBoundary() const
Sampling boundary values instead of cell values.
virtual void print(Ostream &os, int level=0) const
Print information.
virtual const labelList & faceIds() const
Per-face identifier (eg, element Id)
virtual const faceList & faces() const
Faces of surface.
virtual const vectorField & Cf() const
Face centres.
virtual const labelList & zoneIds() const
Per-face zone/region information.
virtual bool hasFaceIds() const
If element ids/order of the original surface are kept.
virtual const scalarField & magSf() const
Face area magnitudes.
virtual bool expire()
Mark the surface as needing an update.
virtual bool needsUpdate() const
Does the surface need an update?
virtual bool update()
Update the surface as required.
virtual const vectorField & Sf() const
Face area vectors.
samplingSource
Types of sampling regions.
@ boundaryFaces
Use nearest boundary face values.
@ insideCells
Surface face within a cell, or trim.
virtual ~sampledMeshedSurface()=default
Destructor.
TypeName("sampledMeshedSurface")
Declare type-name, virtual type (with debug switch)
An abstract class for surfaces with sampling.
const word & name() const noexcept
Name of surface.
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
bool interpolate() const noexcept
Same as isPointData()
A class for managing temporary objects.
Definition: tmp.H:65
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
A class for handling words, derived from Foam::string.
Definition: word.H:68
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
MeshedSurface< face > meshedSurface
Field< vector > vectorField
Specialisation of Field<T> for vector.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73