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 -------------------------------------------------------------------------------
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 
30 #include "dictionary.H"
31 #include "polyMesh.H"
32 #include "polyPatch.H"
33 #include "volFields.H"
34 
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(sampledPatchInternalField, 0);
43  (
44  sampledSurface,
45  sampledPatchInternalField,
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 {
65  mappedPatchBase::offsetModeNames_.getOrDefault
66  (
67  "offsetMode", dict, mappedPatchBase::NORMAL
68  );
69 
70  switch (mode)
71  {
72  case mappedPatchBase::NORMAL:
73  {
74  const scalar distance(dict.get<scalar>("distance"));
75  forAll(patchIDs(), i)
76  {
77  mappers_.set
78  (
79  i,
80  new mappedPatchBase
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 
93  case mappedPatchBase::UNIFORM:
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 
114  case mappedPatchBase::NONUNIFORM:
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 (
142  const interpolation<scalar>& sampler
143 ) const
144 {
145  return sampleOnFaces(sampler);
146 }
147 
148 
150 (
151  const interpolation<vector>& sampler
152 ) const
153 {
154  return sampleOnFaces(sampler);
155 }
156 
157 
159 (
160  const interpolation<sphericalTensor>& sampler
161 ) const
162 {
163  return sampleOnFaces(sampler);
164 }
165 
166 
168 (
169  const interpolation<symmTensor>& sampler
170 ) const
171 {
172  return sampleOnFaces(sampler);
173 }
174 
175 
177 (
178  const interpolation<tensor>& sampler
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 (
206  const interpolation<sphericalTensor>& interpolator
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 // ************************************************************************* //
volFields.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
polyPatch.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::sampledPatch::patchNames
const wordRes & patchNames() const
The selection (word/regex) of patches.
Definition: sampledPatch.H:159
Foam::sampledPatchInternalField::print
virtual void print(Ostream &os, int level=0) const
Print information.
Definition: sampledPatchInternalField.C:231
Foam::sampledPatchInternalField::sample
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample boundary of volume field onto surface faces.
Definition: sampledPatchInternalField.C:141
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
polyMesh.H
Foam::sampledPatch::points
virtual const pointField & points() const
Points of surface.
Definition: sampledPatch.H:232
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::mode
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: MSwindows.C:564
Foam::sampledSurface::interpolate
bool interpolate() const noexcept
Same as isPointData()
Definition: sampledSurface.H:598
Foam::Field< vector >
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
Foam::interpolation< scalar >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::sampledPatchInternalField::sampledPatchInternalField
sampledPatchInternalField(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: sampledPatchInternalField.C:55
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::sampledSurface::name
const word & name() const noexcept
Name of surface.
Definition: sampledSurface.H:322
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
sampledPatchInternalField.H
Foam::Vector< scalar >
dictionary.H
Foam::mappedPatchBase::offsetMode
offsetMode
How to project face centres.
Definition: mappedPatchBase.H:130
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::sampledPatch
A sampledSurface on patches. Non-triangulated by default.
Definition: sampledPatch.H:95
Foam::sampledPatch::faces
virtual const faceList & faces() const
Faces of surface.
Definition: sampledPatch.H:238