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-2018 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"
27 #include "polyMesh.H"
28 #include "demandDrivenData.H"
29 
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(sampledSurface, 0);
36  defineRunTimeSelectionTable(sampledSurface, word);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 void Foam::sampledSurface::makeSf() const
43 {
44  // It is an error to recalculate if the pointer is already set
45  if (SfPtr_)
46  {
48  << "face area vectors already exist"
49  << abort(FatalError);
50  }
51 
52  const faceList& theFaces = faces();
53  SfPtr_ = new vectorField(theFaces.size());
54 
55  vectorField& values = *SfPtr_;
56  forAll(theFaces, facei)
57  {
58  values[facei] = theFaces[facei].area(points());
59  }
60 }
61 
62 
63 void Foam::sampledSurface::makeMagSf() const
64 {
65  // It is an error to recalculate if the pointer is already set
66  if (magSfPtr_)
67  {
69  << "mag face areas already exist"
70  << abort(FatalError);
71  }
72 
73  const faceList& theFaces = faces();
74  magSfPtr_ = new scalarField(theFaces.size());
75 
76  scalarField& values = *magSfPtr_;
77  forAll(theFaces, facei)
78  {
79  values[facei] = theFaces[facei].mag(points());
80  }
81 }
82 
83 
84 void Foam::sampledSurface::makeCf() const
85 {
86  // It is an error to recalculate if the pointer is already set
87  if (CfPtr_)
88  {
90  << "face centres already exist"
91  << abort(FatalError);
92  }
93 
94  const faceList& theFaces = faces();
95  CfPtr_ = new vectorField(theFaces.size());
96 
97  vectorField& values = *CfPtr_;
98  forAll(theFaces, facei)
99  {
100  values[facei] = theFaces[facei].centre(points());
101  }
102 }
103 
104 
105 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
106 
108 {
109  deleteDemandDrivenData(SfPtr_);
110  deleteDemandDrivenData(magSfPtr_);
111  deleteDemandDrivenData(CfPtr_);
112  area_ = -1;
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
117 
119 (
120  const word& name,
121  const polyMesh& mesh,
122  const dictionary& dict
123 )
124 {
125  const word sampleType(dict.lookup("type"));
126 
127  if (debug)
128  {
129  Info<< "Selecting sampledType " << sampleType << endl;
130  }
131 
132  wordConstructorTable::iterator cstrIter =
133  wordConstructorTablePtr_->find(sampleType);
134 
135  if (cstrIter == wordConstructorTablePtr_->end())
136  {
138  << "Unknown sample type "
139  << sampleType << nl << nl
140  << "Valid sample types : " << endl
141  << wordConstructorTablePtr_->sortedToc()
142  << exit(FatalError);
143  }
144 
145  return autoPtr<sampledSurface>(cstrIter()(name, mesh, dict));
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
150 
152 (
153  const word& name,
154  const polyMesh& mesh,
155  const bool interpolate
156 )
157 :
158  name_(name),
159  mesh_(mesh),
160  interpolate_(interpolate),
161  SfPtr_(nullptr),
162  magSfPtr_(nullptr),
163  CfPtr_(nullptr),
164  area_(-1)
165 {}
166 
167 
169 (
170  const word& name,
171  const polyMesh& mesh,
172  const dictionary& dict
173 )
174 :
175  name_(name),
176  mesh_(mesh),
177  interpolate_(dict.lookupOrDefault("interpolate", false)),
178  SfPtr_(nullptr),
179  magSfPtr_(nullptr),
180  CfPtr_(nullptr),
181  area_(-1)
182 {
183  dict.readIfPresent("name", name_);
184 }
185 
186 
187 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
188 
190 {
191  clearGeom();
192 }
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
198 {
199  if (!SfPtr_)
200  {
201  makeSf();
202  }
203 
204  return *SfPtr_;
205 }
206 
207 
209 {
210  if (!magSfPtr_)
211  {
212  makeMagSf();
213  }
214 
215  return *magSfPtr_;
216 }
217 
218 
220 {
221  if (!CfPtr_)
222  {
223  makeCf();
224  }
225 
226  return *CfPtr_;
227 }
228 
229 
230 Foam::scalar Foam::sampledSurface::area() const
231 {
232  if (area_ < 0)
233  {
234  area_ = sum(magSf());
235  reduce(area_, sumOp<scalar>());
236  }
237 
238  return area_;
239 }
240 
241 
243 (
244  const surfaceScalarField& sField
245 ) const
246 {
248  return tmp<scalarField>(nullptr);
249 }
250 
251 
253 (
254  const surfaceVectorField& sField
255 ) const
256 {
258  return tmp<vectorField>(nullptr);
259 }
260 
261 
263 (
264  const surfaceSphericalTensorField& sField
265 ) const
266 {
268  return tmp<sphericalTensorField>(nullptr);
269 }
270 
271 
273 (
274  const surfaceSymmTensorField& sField
275 ) const
276 {
278  return tmp<symmTensorField>(nullptr);
279 }
280 
281 
283 (
284  const surfaceTensorField& sField
285 ) const
286 {
288  return tmp<tensorField>(nullptr);
289 }
290 
291 
293 Foam::sampledSurface::project(const Field<scalar>& field) const
294 {
295  tmp<Field<scalar>> tRes(new Field<scalar>(faces().size()));
296  Field<scalar>& res = tRes.ref();
297 
298  forAll(faces(), facei)
299  {
300  res[facei] = field[facei];
301  }
302 
303  return tRes;
304 }
305 
306 
308 Foam::sampledSurface::project(const Field<vector>& field) const
309 {
310  tmp<Field<scalar>> tRes(new Field<scalar>(faces().size()));
311  project(tRes.ref(), field);
312  return tRes;
313 }
314 
315 
317 Foam::sampledSurface::project(const Field<sphericalTensor>& field) const
318 {
319  tmp<Field<vector>> tRes(new Field<vector>(faces().size()));
320  project(tRes.ref(), field);
321  return tRes;
322 }
323 
324 
326 Foam::sampledSurface::project(const Field<symmTensor>& field) const
327 {
328  tmp<Field<vector>> tRes(new Field<vector>(faces().size()));
329  project(tRes.ref(), field);
330  return tRes;
331 }
332 
333 
335 Foam::sampledSurface::project(const Field<tensor>& field) const
336 {
337  tmp<Field<vector>> tRes(new Field<vector>(faces().size()));
338  project(tRes.ref(), field);
339  return tRes;
340 }
341 
342 
344 {
345  os << type();
346 }
347 
348 
349 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
350 
352 {
353  s.print(os);
354  os.check("Ostream& operator<<(Ostream&, const sampledSurface&");
355  return os;
356 }
357 
358 
359 // ************************************************************************* //
dictionary dict
virtual tmp< scalarField > sample(const volScalarField &) const =0
Sample field on surface.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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:174
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
virtual void print(Ostream &) const
Write.
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:53
static const char nl
Definition: Ostream.H:265
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)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Template functions to aid in the implementation of demand driven data.
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 &)
messageStream Info
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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:366
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