sampledIsoSurfaceCell.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 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 "sampledIsoSurfaceCell.H"
27 #include "dictionary.H"
28 #include "volFields.H"
29 #include "volPointInterpolation.H"
31 #include "fvMesh.H"
32 #include "isoSurfaceCell.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(sampledIsoSurfaceCell, 0);
40  (
41  sampledSurface,
42  sampledIsoSurfaceCell,
43  word,
44  isoSurfaceCell
45  );
46 }
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 bool Foam::sampledIsoSurfaceCell::updateGeometry() const
51 {
52  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
53 
54  // no update needed
55  if (fvm.time().timeIndex() == prevTimeIndex_)
56  {
57  return false;
58  }
59 
60  prevTimeIndex_ = fvm.time().timeIndex();
61 
62  // Clear any stored topo
63  facesPtr_.clear();
64 
65  // Clear derived data
67 
68  // Optionally read volScalarField
69  autoPtr<volScalarField> readFieldPtr_;
70 
71  // 1. see if field in database
72  // 2. see if field can be read
73  const volScalarField* cellFldPtr = NULL;
74  if (fvm.foundObject<volScalarField>(isoField_))
75  {
76  if (debug)
77  {
78  Info<< "sampledIsoSurfaceCell::updateGeometry() : lookup "
79  << isoField_ << endl;
80  }
81 
82  cellFldPtr = &fvm.lookupObject<volScalarField>(isoField_);
83  }
84  else
85  {
86  // Bit of a hack. Read field and store.
87 
88  if (debug)
89  {
90  Info<< "sampledIsoSurfaceCell::updateGeometry() : reading "
91  << isoField_ << " from time " <<fvm.time().timeName()
92  << endl;
93  }
94 
95  readFieldPtr_.reset
96  (
97  new volScalarField
98  (
99  IOobject
100  (
101  isoField_,
102  fvm.time().timeName(),
103  fvm,
106  false
107  ),
108  fvm
109  )
110  );
111 
112  cellFldPtr = readFieldPtr_.operator->();
113  }
114  const volScalarField& cellFld = *cellFldPtr;
115 
116  tmp<pointScalarField> pointFld
117  (
119  );
120 
121  if (average_)
122  {
123  //- From point field and interpolated cell.
124  scalarField cellAvg(fvm.nCells(), scalar(0.0));
125  labelField nPointCells(fvm.nCells(), 0);
126  {
127  for (label pointI = 0; pointI < fvm.nPoints(); pointI++)
128  {
129  const labelList& pCells = fvm.pointCells(pointI);
130 
131  forAll(pCells, i)
132  {
133  label cellI = pCells[i];
134 
135  cellAvg[cellI] += pointFld().internalField()[pointI];
136  nPointCells[cellI]++;
137  }
138  }
139  }
140  forAll(cellAvg, cellI)
141  {
142  cellAvg[cellI] /= nPointCells[cellI];
143  }
144 
145  const isoSurfaceCell iso
146  (
147  fvm,
148  cellAvg,
149  pointFld().internalField(),
150  isoVal_,
151  regularise_
152  );
153 
154  const_cast<sampledIsoSurfaceCell&>
155  (
156  *this
157  ).triSurface::operator=(iso);
158  meshCells_ = iso.meshCells();
159  }
160  else
161  {
162  //- Direct from cell field and point field. Gives bad continuity.
163  const isoSurfaceCell iso
164  (
165  fvm,
166  cellFld.internalField(),
167  pointFld().internalField(),
168  isoVal_,
169  regularise_
170  );
171 
172  const_cast<sampledIsoSurfaceCell&>
173  (
174  *this
175  ).triSurface::operator=(iso);
176  meshCells_ = iso.meshCells();
177  }
178 
179 
180  if (debug)
181  {
182  Pout<< "sampledIsoSurfaceCell::updateGeometry() : constructed iso:"
183  << nl
184  << " regularise : " << regularise_ << nl
185  << " average : " << average_ << nl
186  << " isoField : " << isoField_ << nl
187  << " isoValue : " << isoVal_ << nl
188  << " points : " << points().size() << nl
189  << " tris : " << triSurface::size() << nl
190  << " cut cells : " << meshCells_.size() << endl;
191  }
192 
193  return true;
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
198 
200 (
201  const word& name,
202  const polyMesh& mesh,
203  const dictionary& dict
204 )
205 :
206  sampledSurface(name, mesh, dict),
207  isoField_(dict.lookup("isoField")),
208  isoVal_(readScalar(dict.lookup("isoValue"))),
209  regularise_(dict.lookupOrDefault("regularise", true)),
210  average_(dict.lookupOrDefault("average", true)),
211  zoneKey_(keyType::null),
212  facesPtr_(NULL),
213  prevTimeIndex_(-1),
214  meshCells_(0)
215 {
216 // dict.readIfPresent("zone", zoneKey_);
217 //
218 // if (debug && zoneKey_.size() && mesh.cellZones().findZoneID(zoneKey_) < 0)
219 // {
220 // Info<< "cellZone " << zoneKey_
221 // << " not found - using entire mesh" << endl;
222 // }
223 }
224 
225 
226 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
227 
229 {}
230 
231 
232 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
233 
235 {
236  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
237 
238  return fvm.time().timeIndex() != prevTimeIndex_;
239 }
240 
241 
243 {
244  facesPtr_.clear();
245 
246  // Clear derived data
248 
249  // already marked as expired
250  if (prevTimeIndex_ == -1)
251  {
252  return false;
253  }
254 
255  // force update
256  prevTimeIndex_ = -1;
257  return true;
258 }
259 
260 
262 {
263  return updateGeometry();
264 }
265 
266 
269 (
270  const volScalarField& vField
271 ) const
272 {
273  return sampleField(vField);
274 }
275 
276 
279 (
280  const volVectorField& vField
281 ) const
282 {
283  return sampleField(vField);
284 }
285 
286 
289 (
290  const volSphericalTensorField& vField
291 ) const
292 {
293  return sampleField(vField);
294 }
295 
296 
299 (
300  const volSymmTensorField& vField
301 ) const
302 {
303  return sampleField(vField);
304 }
305 
306 
309 (
310  const volTensorField& vField
311 ) const
312 {
313  return sampleField(vField);
314 }
315 
316 
319 (
320  const interpolation<scalar>& interpolator
321 ) const
322 {
323  return interpolateField(interpolator);
324 }
325 
326 
329 (
330  const interpolation<vector>& interpolator
331 ) const
332 {
333  return interpolateField(interpolator);
334 }
335 
338 (
339  const interpolation<sphericalTensor>& interpolator
340 ) const
341 {
342  return interpolateField(interpolator);
343 }
344 
345 
348 (
349  const interpolation<symmTensor>& interpolator
350 ) const
351 {
352  return interpolateField(interpolator);
353 }
354 
355 
358 (
359  const interpolation<tensor>& interpolator
360 ) const
361 {
362  return interpolateField(interpolator);
363 }
364 
365 
367 {
368  os << "sampledIsoSurfaceCell: " << name() << " :"
369  << " field:" << isoField_
370  << " value:" << isoVal_;
371  //<< " faces:" << faces().size() // possibly no geom yet
372  //<< " points:" << points().size();
373 }
374 
375 
376 // ************************************************************************* //
virtual bool update()
Update the surface as required.
const pointField & points
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
An abstract class for surfaces with sampling.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const volPointInterpolation & New(const fvMesh &mesh)
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
A class for handling words, derived from string.
Definition: word.H:59
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
friend Ostream & operator(Ostream &, const UList< T > &)
messageStream Info
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
label timeIndex() const
Return current time index.
Definition: TimeState.C:73
static const keyType null
An empty keyType.
Definition: keyType.H:75
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
virtual void clearGeom() const
Macros for easy insertion into run-time selection tables.
sampledIsoSurfaceCell(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
bool interpolate() const
Interpolation requested for surface.
virtual ~sampledIsoSurfaceCell()
Destructor.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
virtual bool expire()
Mark the surface as needing an update.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
List< label > labelList
A List of labels.
Definition: labelList.H:56
label size() const
Return the number of elements in the UList.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
addNamedToRunTimeSelectionTable(GAMGProcAgglomeration, noneGAMGProcAgglomeration, GAMGAgglomeration, none)
A class for managing temporary objects.
Definition: PtrList.H:118
virtual void print(Ostream &) const
Write.
conserve internalField()+
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
virtual bool needsUpdate() const
Does the surface need an update?