Loading...
Searching...
No Matches
sampledPatchInternalField.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-2013 OpenFOAM Foundation
9 Copyright (C) 2018-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
27\*---------------------------------------------------------------------------*/
28
30#include "dictionary.H"
31#include "polyMesh.H"
32#include "polyPatch.H"
33#include "volFields.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
43 (
46 word,
47 patchInternalField
48 );
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const word& name,
57 const polyMesh& mesh,
58 const dictionary& dict
59)
60:
62 mappers_(patchIDs().size())
63{
66 (
67 "offsetMode", dict, mappedPatchBase::NORMAL
68 );
69
70 switch (mode)
71 {
73 {
74 const scalar distance(dict.get<scalar>("distance"));
75 forAll(patchIDs(), i)
76 {
77 mappers_.set
78 (
79 i,
81 (
82 mesh.boundaryMesh()[patchIDs()[i]],
83 mesh.name(), // sampleRegion
84 mappedPatchBase::NEARESTCELL, // sampleMode
85 word::null, // samplePatch
86 -distance // sample inside my domain
87 )
88 );
89 }
90 }
91 break;
92
94 {
95 const point offset(dict.get<point>("offset"));
96 forAll(patchIDs(), i)
97 {
98 mappers_.set
99 (
100 i,
101 new mappedPatchBase
102 (
103 mesh.boundaryMesh()[patchIDs()[i]],
104 mesh.name(), // sampleRegion
105 mappedPatchBase::NEARESTCELL, // sampleMode
106 word::null, // samplePatch
107 offset // sample inside my domain
108 )
109 );
110 }
111 }
112 break;
113
115 {
116 const pointField offsets(dict.get<pointField>("offsets"));
117 forAll(patchIDs(), i)
118 {
119 mappers_.set
120 (
121 i,
122 new mappedPatchBase
123 (
124 mesh.boundaryMesh()[patchIDs()[i]],
125 mesh.name(), // sampleRegion
126 mappedPatchBase::NEARESTCELL, // sampleMode
127 word::null, // samplePatch
128 offsets // sample inside my domain
129 )
130 );
131 }
132 }
133 break;
134 }
135}
136
137
138// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139
141(
143) const
144{
145 return sampleOnFaces(sampler);
146}
147
148
150(
152) const
153{
154 return sampleOnFaces(sampler);
155}
156
157
159(
161) const
162{
163 return sampleOnFaces(sampler);
164}
165
166
168(
170) const
171{
172 return sampleOnFaces(sampler);
173}
174
175
177(
179) const
180{
181 return sampleOnFaces(sampler);
182}
183
184
186(
187 const interpolation<scalar>& interpolator
188) const
189{
190 return sampleOnPoints(interpolator);
191}
192
193
195(
196 const interpolation<vector>& interpolator
197) const
198{
199 return sampleOnPoints(interpolator);
200}
201
202
205(
207) const
208{
209 return sampleOnPoints(interpolator);
210}
211
212
214(
215 const interpolation<symmTensor>& interpolator
216) const
217{
218 return sampleOnPoints(interpolator);
219}
220
221
223(
224 const interpolation<tensor>& interpolator
225) const
226{
227 return sampleOnPoints(interpolator);
228}
229
230
232{
233 os << "sampledPatchInternalField: " << name() << " :"
234 << " patches:" << patchNames();
235
236 if (level)
237 {
238 os << " faces:" << faces().size()
239 << " points:" << points().size();
240 }
241}
242
243
244// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
labelList patchIDs
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base class for volume field interpolation.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
@ NEARESTCELL
nearest cell containing sample
offsetMode
How to project face centres.
@ NORMAL
use face normal + distance
@ UNIFORM
single offset vector
@ NONUNIFORM
per-face offset vector
static const Enum< offsetMode > offsetModeNames_
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Variation of sampledPatch that samples the internalField (at a given normal distance from the patch) ...
sampledPatchInternalField(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
virtual void print(Ostream &os, int level=0) const
Print information.
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample boundary of volume field onto surface faces.
A sampledSurface on patches. Non-triangulated by default.
const labelList & patchIDs() const
The patches selected.
virtual const faceList & faces() const
Faces of surface.
sampledPatch(const word &name, const polyMesh &mesh, const UList< wordRe > &patchNames, const bool triangulate=false)
Construct from components.
An abstract class for surfaces with sampling.
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:75
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
auto & name
const pointField & points
Namespace for OpenFOAM.
scalar distance(const vector &p1, const vector &p2)
Definition curveTools.C:12
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition POSIX.C:775
vector point
Point is a vector.
Definition point.H:37
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
wordList patchNames(nPatches)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299