sampledPatchInternalField.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-2019 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 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace sampledSurfaces
34 {
35 
36  defineTypeNameAndDebug(patchInternalField, 0);
37  addToRunTimeSelectionTable(sampledSurface, patchInternalField, word);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const word& name,
47  const polyMesh& mesh,
48  const dictionary& dict
49 )
50 :
51  patch(name, mesh, dict),
52  mappers_(patchIDs().size())
53 {
55  if (dict.found("offsetMode"))
56  {
58  (
59  dict.lookup("offsetMode")
60  );
61  }
62 
63  switch (mode)
64  {
66  {
67  const scalar distance = dict.lookup<scalar>("distance");
68  forAll(patchIDs(), i)
69  {
70  mappers_.set
71  (
72  i,
73  new mappedPatchBase
74  (
75  mesh.boundaryMesh()[patchIDs()[i]],
76  mesh.name(), // sampleRegion
77  mappedPatchBase::NEARESTCELL, // sampleMode
78  word::null, // samplePatch
79  -distance // sample inside my domain
80  )
81  );
82  }
83  }
84  break;
85 
87  {
88  const point offset(dict.lookup("offset"));
89  forAll(patchIDs(), i)
90  {
91  mappers_.set
92  (
93  i,
94  new mappedPatchBase
95  (
96  mesh.boundaryMesh()[patchIDs()[i]],
97  mesh.name(), // sampleRegion
98  mappedPatchBase::NEARESTCELL, // sampleMode
99  word::null, // samplePatch
100  offset // sample inside my domain
101  )
102  );
103  }
104  }
105  break;
106 
108  {
109  const pointField offsets(dict.lookup("offsets"));
110  forAll(patchIDs(), i)
111  {
112  mappers_.set
113  (
114  i,
115  new mappedPatchBase
116  (
117  mesh.boundaryMesh()[patchIDs()[i]],
118  mesh.name(), // sampleRegion
119  mappedPatchBase::NEARESTCELL, // sampleMode
120  word::null, // samplePatch
121  offsets // sample inside my domain
122  )
123  );
124  }
125  }
126  break;
127  }
128 }
129 
130 
131 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
132 
134 {}
135 
136 
137 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 
141 (
142  const volScalarField& vField
143 ) const
144 {
145  return sampleField(vField);
146 }
147 
148 
151 (
152  const volVectorField& vField
153 ) const
154 {
155  return sampleField(vField);
156 }
157 
160 (
161  const volSphericalTensorField& vField
162 ) const
163 {
164  return sampleField(vField);
165 }
166 
167 
170 (
171  const volSymmTensorField& vField
172 ) const
173 {
174  return sampleField(vField);
175 }
176 
177 
180 (
181  const volTensorField& vField
182 ) const
183 {
184  return sampleField(vField);
185 }
186 
187 
190 (
191  const interpolation<scalar>& interpolator
192 ) const
193 {
194  return interpolateField(interpolator);
195 }
196 
197 
200 (
201  const interpolation<vector>& interpolator
202 ) const
203 {
204  return interpolateField(interpolator);
205 }
206 
207 
210 (
211  const interpolation<sphericalTensor>& interpolator
212 ) const
213 {
214  return interpolateField(interpolator);
215 }
216 
217 
220 (
221  const interpolation<symmTensor>& interpolator
222 ) const
223 {
224  return interpolateField(interpolator);
225 }
226 
227 
230 (
231  const interpolation<tensor>& interpolator
232 ) const
233 {
234  return interpolateField(interpolator);
235 }
236 
237 
239 {
240  os << "patchInternalField: " << name() << " :"
241  << " patches:" << patchNames()
242  << " faces:" << faces().size()
243  << " points:" << points().size();
244 }
245 
246 
247 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
mode_t mode(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file mode.
Definition: POSIX.C:461
const word & name() const
Return name.
Definition: IOobject.H:315
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
patchInternalField(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
defineTypeNameAndDebug(distanceSurface, 0)
A sampledSurface on patches. Non-triangulated by default.
Definition: sampledPatch.H:89
bool interpolate() const
Interpolation requested for surface.
virtual void print(Ostream &) const
Write.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Macros for easy insertion into run-time selection tables.
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
wordList patchNames(nPatches)
static const word null
An empty word.
Definition: word.H:77
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
offsetMode
How to project face centres.
addToRunTimeSelectionTable(sampledSurface, distanceSurface, word)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:61
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
static const NamedEnum< offsetMode, 3 > offsetModeNames_
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864