sampledSurface.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 
26 #include "sampledSurface.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(sampledSurface, 0);
34  defineRunTimeSelectionTable(sampledSurface, word);
35 }
36 
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 void Foam::sampledSurface::makeSf() const
41 {
42  // It is an error to recalculate if the pointer is already set
43  if (SfPtr_)
44  {
46  << "face area vectors already exist"
47  << abort(FatalError);
48  }
49 
50  const faceList& theFaces = faces();
51  SfPtr_ = new vectorField(theFaces.size());
52 
53  vectorField& values = *SfPtr_;
54  forAll(theFaces, facei)
55  {
56  values[facei] = theFaces[facei].area(points());
57  }
58 }
59 
60 
61 void Foam::sampledSurface::makeMagSf() const
62 {
63  // It is an error to recalculate if the pointer is already set
64  if (magSfPtr_)
65  {
67  << "mag face areas already exist"
68  << abort(FatalError);
69  }
70 
71  const faceList& theFaces = faces();
72  magSfPtr_ = new scalarField(theFaces.size());
73 
74  scalarField& values = *magSfPtr_;
75  forAll(theFaces, facei)
76  {
77  values[facei] = theFaces[facei].mag(points());
78  }
79 }
80 
81 
82 void Foam::sampledSurface::makeCf() const
83 {
84  // It is an error to recalculate if the pointer is already set
85  if (CfPtr_)
86  {
88  << "face centres already exist"
89  << abort(FatalError);
90  }
91 
92  const faceList& theFaces = faces();
93  CfPtr_ = new vectorField(theFaces.size());
94 
95  vectorField& values = *CfPtr_;
96  forAll(theFaces, facei)
97  {
98  values[facei] = theFaces[facei].centre(points());
99  }
100 }
101 
102 
103 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
104 
106 {
107  deleteDemandDrivenData(SfPtr_);
108  deleteDemandDrivenData(magSfPtr_);
109  deleteDemandDrivenData(CfPtr_);
110  area_ = -1;
111 }
112 
113 
114 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
115 
117 (
118  const word& name,
119  const polyMesh& mesh,
120  const dictionary& dict
121 )
122 {
123  const word sampleType(dict.lookup("type"));
124 
125  if (debug)
126  {
127  Info<< "Selecting sampledType " << sampleType << endl;
128  }
129 
130  wordConstructorTable::iterator cstrIter =
131  wordConstructorTablePtr_->find(sampleType);
132 
133  if (cstrIter == wordConstructorTablePtr_->end())
134  {
136  << "Unknown sample type "
137  << sampleType << nl << nl
138  << "Valid sample types : " << endl
139  << wordConstructorTablePtr_->sortedToc()
140  << exit(FatalError);
141  }
142 
143  return autoPtr<sampledSurface>(cstrIter()(name, mesh, dict));
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
148 
150 (
151  const word& name,
152  const polyMesh& mesh,
153  const bool interpolate
154 )
155 :
156  name_(name),
157  mesh_(mesh),
158  interpolate_(interpolate),
159  SfPtr_(nullptr),
160  magSfPtr_(nullptr),
161  CfPtr_(nullptr),
162  area_(-1)
163 {}
164 
165 
167 (
168  const word& name,
169  const polyMesh& mesh,
170  const dictionary& dict
171 )
172 :
173  name_(name),
174  mesh_(mesh),
175  interpolate_(dict.lookupOrDefault("interpolate", false)),
176  SfPtr_(nullptr),
177  magSfPtr_(nullptr),
178  CfPtr_(nullptr),
179  area_(-1)
180 {
181  dict.readIfPresent("name", name_);
182 }
183 
184 
185 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
186 
188 {
189  clearGeom();
190 }
191 
192 
193 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
194 
196 {
197  if (!SfPtr_)
198  {
199  makeSf();
200  }
201 
202  return *SfPtr_;
203 }
204 
205 
207 {
208  if (!magSfPtr_)
209  {
210  makeMagSf();
211  }
212 
213  return *magSfPtr_;
214 }
215 
216 
218 {
219  if (!CfPtr_)
220  {
221  makeCf();
222  }
223 
224  return *CfPtr_;
225 }
226 
227 
228 Foam::scalar Foam::sampledSurface::area() const
229 {
230  if (area_ < 0)
231  {
232  area_ = sum(magSf());
233  reduce(area_, sumOp<scalar>());
234  }
235 
236  return area_;
237 }
238 
239 
241 (
242  const surfaceScalarField& sField
243 ) const
244 {
246  return tmp<scalarField>(nullptr);
247 }
248 
249 
251 (
252  const surfaceVectorField& sField
253 ) const
254 {
256  return tmp<vectorField>(nullptr);
257 }
258 
259 
261 (
262  const surfaceSphericalTensorField& sField
263 ) const
264 {
266  return tmp<sphericalTensorField>(nullptr);
267 }
268 
269 
271 (
272  const surfaceSymmTensorField& sField
273 ) const
274 {
276  return tmp<symmTensorField>(nullptr);
277 }
278 
279 
281 (
282  const surfaceTensorField& sField
283 ) const
284 {
286  return tmp<tensorField>(nullptr);
287 }
288 
289 
291 Foam::sampledSurface::project(const Field<scalar>& field) const
292 {
293  tmp<Field<scalar>> tRes(new Field<scalar>(faces().size()));
294  Field<scalar>& res = tRes.ref();
295 
296  forAll(faces(), facei)
297  {
298  res[facei] = field[facei];
299  }
300 
301  return tRes;
302 }
303 
304 
306 Foam::sampledSurface::project(const Field<vector>& field) const
307 {
308  tmp<Field<scalar>> tRes(new Field<scalar>(faces().size()));
309  project(tRes.ref(), field);
310  return tRes;
311 }
312 
313 
315 Foam::sampledSurface::project(const Field<sphericalTensor>& field) const
316 {
317  tmp<Field<vector>> tRes(new Field<vector>(faces().size()));
318  project(tRes.ref(), field);
319  return tRes;
320 }
321 
322 
324 Foam::sampledSurface::project(const Field<symmTensor>& field) const
325 {
326  tmp<Field<vector>> tRes(new Field<vector>(faces().size()));
327  project(tRes.ref(), field);
328  return tRes;
329 }
330 
331 
333 Foam::sampledSurface::project(const Field<tensor>& field) const
334 {
335  tmp<Field<vector>> tRes(new Field<vector>(faces().size()));
336  project(tRes.ref(), field);
337  return tRes;
338 }
339 
340 
342 {
343  os << type();
344 }
345 
346 
347 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
348 
350 {
351  s.print(os);
352  os.check("Ostream& operator<<(Ostream&, const sampledSurface&");
353  return os;
354 }
355 
356 
357 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual ~sampledSurface()
Destructor.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
An abstract class for surfaces with sampling.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void print(Ostream &) const
Write.
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dynamicFvMesh & mesh
virtual void clearGeom() const
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
sampledSurface(const word &name, const polyMesh &, const bool interpolate=false)
Construct from name, mesh.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
scalar area() const
The total surface area.
static autoPtr< sampledSurface > New(const word &name, const polyMesh &, const dictionary &)
Return a reference to the selected surface.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
virtual const vectorField & Sf() const
Return face area vectors.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
tmp< Field< vector > > project(const Field< sphericalTensor > &) const
Project field onto surface.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual const vectorField & Cf() const
Return face centres as vectorField.
virtual const scalarField & magSf() const
Return face area magnitudes.
Ostream & operator<<(Ostream &, const ensightPart &)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
messageStream Info
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
rDeltaTY field()
A class for managing temporary objects.
Definition: PtrList.H:53
void deleteDemandDrivenData(DataPtr &dataPtr)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844
virtual tmp< scalarField > sample(const volScalarField &) const =0
Sample field on surface.