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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "sampledIsoSurfaceCell.H"
30 #include "dictionary.H"
31 #include "volFields.H"
32 #include "volPointInterpolation.H"
34 #include "fvMesh.H"
35 #include "isoSurfaceCell.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(sampledIsoSurfaceCell, 0);
43  (
44  sampledSurface,
45  sampledIsoSurfaceCell,
46  word,
47  isoSurfaceCell
48  );
49 }
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 bool Foam::sampledIsoSurfaceCell::updateGeometry() const
54 {
55  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
56 
57  // No update needed
58  if (fvm.time().timeIndex() == prevTimeIndex_)
59  {
60  return false;
61  }
62 
63  prevTimeIndex_ = fvm.time().timeIndex();
64 
65  // Clear derived data
67 
68  // Use field from database, or try to read it in
69 
70  const auto* cellFldPtr = fvm.findObject<volScalarField>(isoField_);
71 
72  if (debug)
73  {
74  if (cellFldPtr)
75  {
76  InfoInFunction << "Lookup " << isoField_ << endl;
77  }
78  else
79  {
81  << "Reading " << isoField_
82  << " from time " << fvm.time().timeName()
83  << endl;
84  }
85  }
86 
87  // For holding the volScalarField read in.
88  autoPtr<volScalarField> fieldReadPtr;
89 
90  if (!cellFldPtr)
91  {
92  // Bit of a hack. Read field and store.
93 
94  fieldReadPtr = autoPtr<volScalarField>::New
95  (
96  IOobject
97  (
98  isoField_,
99  fvm.time().timeName(),
100  fvm,
103  false
104  ),
105  fvm
106  );
107  }
108 
109  const volScalarField& cellFld =
110  (fieldReadPtr.valid() ? *fieldReadPtr : *cellFldPtr);
111 
112  auto tpointFld = volPointInterpolation::New(fvm).interpolate(cellFld);
113 
114  if (average_)
115  {
116  //- From point field and interpolated cell.
117  scalarField cellAvg(fvm.nCells(), Zero);
118  labelField nPointCells(fvm.nCells(), Zero);
119 
120  for (label pointi = 0; pointi < fvm.nPoints(); ++pointi)
121  {
122  const scalar& val = tpointFld().primitiveField()[pointi];
123  const labelList& pCells = fvm.pointCells(pointi);
124 
125  for (const label celli : pCells)
126  {
127  cellAvg[celli] += val;
128  ++nPointCells[celli];
129  }
130  }
131  forAll(cellAvg, celli)
132  {
133  cellAvg[celli] /= nPointCells[celli];
134  }
135 
136  isoSurfaceCell surf
137  (
138  fvm,
139  cellAvg,
140  tpointFld().primitiveField(),
141  isoVal_,
142  filter_,
143  bounds_
144  );
145 
146  const_cast<sampledIsoSurfaceCell&>
147  (
148  *this
149  ).transfer(static_cast<meshedSurface&>(surf));
150  meshCells_.transfer(surf.meshCells());
151  }
152  else
153  {
154  //- Direct from cell field and point field. Gives bad continuity.
155  isoSurfaceCell surf
156  (
157  fvm,
158  cellFld.primitiveField(),
159  tpointFld().primitiveField(),
160  isoVal_,
161  filter_,
162  bounds_
163  );
164 
165  const_cast<sampledIsoSurfaceCell&>
166  (
167  *this
168  ).transfer(static_cast<meshedSurface&>(surf));
169  meshCells_.transfer(surf.meshCells());
170  }
171 
172 
173  if (debug)
174  {
175  Pout<< "sampledIsoSurfaceCell::updateGeometry() : constructed iso:"
176  << nl
177  << " filter : " << Switch(bool(filter_)) << nl
178  << " average : " << Switch(average_) << nl
179  << " isoField : " << isoField_ << nl
180  << " isoValue : " << isoVal_ << nl
181  << " bounds : " << bounds_ << nl
182  << " points : " << points().size() << nl
183  << " faces : " << MeshStorage::size() << nl
184  << " cut cells : " << meshCells_.size() << endl;
185  }
186 
187  return true;
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
192 
194 (
195  const word& name,
196  const polyMesh& mesh,
197  const dictionary& dict
198 )
199 :
201  MeshStorage(),
202  isoField_(dict.get<word>("isoField")),
203  isoVal_(dict.get<scalar>("isoValue")),
204  filter_
205  (
207  (
208  dict,
210  )
211  ),
212  average_(dict.getOrDefault("average", true)),
213  bounds_(dict.getOrDefault("bounds", boundBox::invertedBox)),
214  prevTimeIndex_(-1),
215  meshCells_()
216 {}
217 
218 
219 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
220 
222 {}
223 
224 
225 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
226 
228 {
229  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
230 
231  return fvm.time().timeIndex() != prevTimeIndex_;
232 }
233 
234 
236 {
237  // Clear derived data
239 
240  // Already marked as expired
241  if (prevTimeIndex_ == -1)
242  {
243  return false;
244  }
245 
246  // Force update
247  prevTimeIndex_ = -1;
248  return true;
249 }
250 
251 
253 {
254  return updateGeometry();
255 }
256 
257 
260 (
261  const interpolation<scalar>& sampler
262 ) const
263 {
264  return sampleOnFaces(sampler);
265 }
266 
267 
270 (
271  const interpolation<vector>& sampler
272 ) const
273 {
274  return sampleOnFaces(sampler);
275 }
276 
277 
280 (
281  const interpolation<sphericalTensor>& sampler
282 ) const
283 {
284  return sampleOnFaces(sampler);
285 }
286 
287 
290 (
291  const interpolation<symmTensor>& sampler
292 ) const
293 {
294  return sampleOnFaces(sampler);
295 }
296 
297 
300 (
301  const interpolation<tensor>& sampler
302 ) const
303 {
304  return sampleOnFaces(sampler);
305 }
306 
307 
310 (
311  const interpolation<scalar>& interpolator
312 ) const
313 {
314  return sampleOnPoints(interpolator);
315 }
316 
317 
320 (
321  const interpolation<vector>& interpolator
322 ) const
323 {
324  return sampleOnPoints(interpolator);
325 }
326 
329 (
330  const interpolation<sphericalTensor>& interpolator
331 ) const
332 {
333  return sampleOnPoints(interpolator);
334 }
335 
336 
339 (
340  const interpolation<symmTensor>& interpolator
341 ) const
342 {
343  return sampleOnPoints(interpolator);
344 }
345 
346 
349 (
350  const interpolation<tensor>& interpolator
351 ) const
352 {
353  return sampleOnPoints(interpolator);
354 }
355 
356 
358 {
359  os << "sampledIsoSurfaceCell: " << name() << " :"
360  << " field:" << isoField_
361  << " value:" << isoVal_;
362  //<< " faces:" << faces().size() // possibly no geom yet
363  //<< " points:" << points().size();
364 }
365 
366 
367 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::sampledIsoSurfaceCell::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledIsoSurfaceCell.C:235
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::sampledIsoSurfaceCell::print
virtual void print(Ostream &) const
Write.
Definition: sampledIsoSurfaceCell.C:357
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::val
label ListType::const_reference val
Definition: ListOps.H:407
Foam::isoSurfaceBase::filterType::DIAGCELL
Remove points from pyramid edges and face-diagonals.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
isoSurfaceCell.H
Foam::boundBox::invertedBox
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Foam::sampledIsoSurfaceCell::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledIsoSurfaceCell.C:227
Foam::MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation >::New
static const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::sampledIsoSurfaceCell::sample
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
Definition: sampledIsoSurfaceCell.C:260
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::volPointInterpolation::interpolate
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::MeshedSurface< face >::transfer
void transfer(pointField &pointLst, List< face > &faceLst)
Transfer the components.
Definition: MeshedSurface.C:1128
Foam::label
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:62
Foam::labelField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:436
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
Foam::sampledIsoSurfaceCell::points
virtual const pointField & points() const
Points of surface.
Definition: sampledIsoSurfaceCell.H:219
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:120
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::interpolation< scalar >
Foam::sampledIsoSurfaceCell::update
virtual bool update()
Update the surface as required.
Definition: sampledIsoSurfaceCell.C:252
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
volPointInterpolation.H
sampledIsoSurfaceCell.H
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::meshedSurface
MeshedSurface< face > meshedSurface
Definition: MeshedSurfacesFwd.H:41
Foam::sampledIsoSurfaceCell::~sampledIsoSurfaceCell
virtual ~sampledIsoSurfaceCell()
Destructor.
Definition: sampledIsoSurfaceCell.C:221
dictionary.H
Foam::isoSurfaceBase::getFilterType
static filterType getFilterType(const dictionary &dict, const isoSurfaceBase::filterType deflt)
Get 'regularise' as bool or enumeration.
Definition: isoSurfaceBase.C:62
Foam::sampledSurface::clearGeom
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
Definition: sampledSurface.C:55
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::sampledSurface::mesh
const polyMesh & mesh() const
Access to the underlying mesh.
Definition: sampledSurface.H:302
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::TimeState::timeIndex
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:36
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::MeshedSurface< face >
Foam::MeshedSurface< face >::size
label size() const
The surface size is the number of faces.
Definition: MeshedSurface.H:356
Foam::sampledSurface::interpolate
bool interpolate() const
Interpolation to nodes requested for surface.
Definition: sampledSurface.H:326
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::sampledIsoSurfaceCell::sampledIsoSurfaceCell
sampledIsoSurfaceCell(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: sampledIsoSurfaceCell.C:194
Foam::IOobject::MUST_READ
Definition: IOobject.H:120