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 | Copyright (C) 2011-2013 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 
27 #include "dictionary.H"
28 #include "polyMesh.H"
29 #include "polyPatch.H"
30 #include "volFields.H"
31 
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(sampledPatchInternalField, 0);
40  (
41  sampledSurface,
42  sampledPatchInternalField,
43  word,
44  patchInternalField
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& name,
54  const polyMesh& mesh,
55  const dictionary& dict
56 )
57 :
58  sampledPatch(name, mesh, dict),
59  mappers_(patchIDs().size())
60 {
62  if (dict.found("offsetMode"))
63  {
65  (
66  dict.lookup("offsetMode")
67  );
68  }
69 
70  switch (mode)
71  {
73  {
74  const scalar distance = readScalar(dict.lookup("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 
94  {
95  const point offset(dict.lookup("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.lookup("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 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
139 
141 {}
142 
143 
144 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145 
147 (
148  const volScalarField& vField
149 ) const
150 {
151  return sampleField(vField);
152 }
153 
154 
156 (
157  const volVectorField& vField
158 ) const
159 {
160  return sampleField(vField);
161 }
162 
164 (
165  const volSphericalTensorField& vField
166 ) const
167 {
168  return sampleField(vField);
169 }
170 
171 
173 (
174  const volSymmTensorField& vField
175 ) const
176 {
177  return sampleField(vField);
178 }
179 
180 
182 (
183  const volTensorField& vField
184 ) const
185 {
186  return sampleField(vField);
187 }
188 
189 
191 (
192  const interpolation<scalar>& interpolator
193 ) const
194 {
195  return interpolateField(interpolator);
196 }
197 
198 
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 
219 (
220  const interpolation<symmTensor>& interpolator
221 ) const
222 {
223  return interpolateField(interpolator);
224 }
225 
226 
228 (
229  const interpolation<tensor>& interpolator
230 ) const
231 {
232  return interpolateField(interpolator);
233 }
234 
235 
237 {
238  os << "sampledPatchInternalField: " << name() << " :"
239  << " patches:" << patchNames()
240  << " faces:" << faces().size()
241  << " points:" << points().size();
242 }
243 
244 
245 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const word & name() const
Return name.
Definition: IOobject.H:291
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
bool interpolate() const
Interpolation requested for surface.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Macros for easy insertion into run-time selection tables.
A sampledSurface on patches. Non-triangulated by default.
Definition: sampledPatch.H:49
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
mode_t mode(const fileName &, const bool followLink=true)
Return the file mode.
Definition: POSIX.C:467
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 ...
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
defineTypeNameAndDebug(combustionModel, 0)
offsetMode
How to project face centres.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
addNamedToRunTimeSelectionTable(GAMGProcAgglomeration, noneGAMGProcAgglomeration, GAMGAgglomeration, none)
virtual ~sampledPatchInternalField()
Destructor.
virtual void print(Ostream &) const
Write.
sampledPatchInternalField(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
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:74
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:576