sampledIsoSurfaceTopo.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) 2018 OpenFOAM Foundation
9  Copyright (C) 2018-2019 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 "sampledIsoSurfaceTopo.H"
30 #include "dictionary.H"
31 #include "volFields.H"
32 #include "volPointInterpolation.H"
34 #include "fvMesh.H"
35 #include "isoSurfaceTopo.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(sampledIsoSurfaceTopo, 0);
43  (
44  sampledSurface,
45  sampledIsoSurfaceTopo,
46  word,
47  isoSurfaceTopo
48  );
49 }
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 bool Foam::sampledIsoSurfaceTopo::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  //- Direct from cell field and point field. Gives bad continuity.
115  isoSurfaceTopo surf
116  (
117  fvm,
118  cellFld.primitiveField(),
119  tpointFld().primitiveField(),
120  isoVal_,
121  filter_
122  );
123 
124  MeshedSurface<face>& mySurface = const_cast<sampledIsoSurfaceTopo&>(*this);
125 
126  mySurface.transfer(static_cast<meshedSurface&>(surf));
127  meshCells_ = std::move(surf.meshCells());
128 
129  // triangulate uses remapFaces()
130  // - this is somewhat less efficient since it recopies the faces
131  // that we just created, but we probably don't want to do this
132  // too often anyhow.
133  if (triangulate_)
134  {
136  mySurface.triangulate(faceMap);
137  meshCells_ = UIndirectList<label>(meshCells_, faceMap)();
138  }
139 
140  if (debug)
141  {
142  Pout<< "sampledIsoSurfaceTopo::updateGeometry() : constructed iso:"
143  << nl
144  << " filter : " << isoSurfaceBase::filterNames[filter_]
145  << nl
146  << " triangulate : " << Switch(triangulate_) << nl
147  << " isoField : " << isoField_ << nl
148  << " isoValue : " << isoVal_ << nl
149  << " points : " << points().size() << nl
150  << " faces : " << MeshStorage::size() << nl
151  << " cut cells : " << meshCells_.size() << endl;
152  }
153 
154  return true;
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
159 
161 (
162  const word& name,
163  const polyMesh& mesh,
164  const dictionary& dict
165 )
166 :
168  MeshStorage(),
169  isoField_(dict.get<word>("isoField")),
170  isoVal_(dict.get<scalar>("isoValue")),
171  filter_
172  (
174  (
175  dict,
177  )
178  ),
179  triangulate_(dict.getOrDefault("triangulate", false)),
180  prevTimeIndex_(-1),
181  meshCells_()
182 {
183  if (triangulate_ && filter_ == isoSurfaceBase::filterType::NONE)
184  {
186  << "Cannot triangulate without a regularise filter" << nl
187  << exit(FatalIOError);
188  }
189 }
190 
191 
192 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
193 
195 {}
196 
197 
198 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
199 
201 {
202  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
203 
204  return fvm.time().timeIndex() != prevTimeIndex_;
205 }
206 
207 
209 {
210  // Clear derived data
212 
213  // Already marked as expired
214  if (prevTimeIndex_ == -1)
215  {
216  return false;
217  }
218 
219  // Force update
220  prevTimeIndex_ = -1;
221  return true;
222 }
223 
224 
226 {
227  return updateGeometry();
228 }
229 
230 
233 (
234  const interpolation<scalar>& sampler
235 ) const
236 {
237  return sampleOnFaces(sampler);
238 }
239 
240 
243 (
244  const interpolation<vector>& sampler
245 ) const
246 {
247  return sampleOnFaces(sampler);
248 }
249 
250 
253 (
254  const interpolation<sphericalTensor>& sampler
255 ) const
256 {
257  return sampleOnFaces(sampler);
258 }
259 
260 
263 (
264  const interpolation<symmTensor>& sampler
265 ) const
266 {
267  return sampleOnFaces(sampler);
268 }
269 
270 
273 (
274  const interpolation<tensor>& sampler
275 ) const
276 {
277  return sampleOnFaces(sampler);
278 }
279 
280 
283 (
284  const interpolation<scalar>& interpolator
285 ) const
286 {
287  return sampleOnPoints(interpolator);
288 }
289 
290 
293 (
294  const interpolation<vector>& interpolator
295 ) const
296 {
297  return sampleOnPoints(interpolator);
298 }
299 
302 (
303  const interpolation<sphericalTensor>& interpolator
304 ) const
305 {
306  return sampleOnPoints(interpolator);
307 }
308 
309 
312 (
313  const interpolation<symmTensor>& interpolator
314 ) const
315 {
316  return sampleOnPoints(interpolator);
317 }
318 
319 
322 (
323  const interpolation<tensor>& interpolator
324 ) const
325 {
326  return sampleOnPoints(interpolator);
327 }
328 
329 
331 {
332  os << "sampledIsoSurfaceTopo: " << name() << " :"
333  << " field:" << isoField_
334  << " value:" << isoVal_;
335  //<< " faces:" << faces().size() // possibly no geom yet
336  //<< " points:" << points().size();
337 }
338 
339 
340 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
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::isoSurfaceBase::filterType::NONE
No filtering.
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:94
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
isoSurfaceTopo.H
Foam::sampledIsoSurfaceTopo::update
virtual bool update()
Update the surface as required.
Definition: sampledIsoSurfaceTopo.C:225
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
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::FatalIOError
IOerror FatalIOError
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::isoSurfaceBase::filterNames
static const Enum< filterType > filterNames
Names for the filtering types.
Definition: isoSurfaceBase.H:103
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
Foam::sampledIsoSurfaceTopo::print
virtual void print(Ostream &) const
Write.
Definition: sampledIsoSurfaceTopo.C:330
Foam::sampledIsoSurfaceTopo::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledIsoSurfaceTopo.C:200
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
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 >
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::sampledIsoSurfaceTopo::sampledIsoSurfaceTopo
sampledIsoSurfaceTopo(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: sampledIsoSurfaceTopo.C:161
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::sampledIsoSurfaceTopo::points
virtual const pointField & points() const
Points of surface.
Definition: sampledIsoSurfaceTopo.H:205
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
volPointInterpolation.H
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::sampledIsoSurfaceTopo::sample
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
Definition: sampledIsoSurfaceTopo.C:233
Foam::meshedSurface
MeshedSurface< face > meshedSurface
Definition: MeshedSurfacesFwd.H:41
Foam::sampledIsoSurfaceTopo::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledIsoSurfaceTopo.C:208
sampledIsoSurfaceTopo.H
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
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
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::sampledIsoSurfaceTopo::~sampledIsoSurfaceTopo
virtual ~sampledIsoSurfaceTopo()
Destructor.
Definition: sampledIsoSurfaceTopo.C:194
Foam::IOobject::MUST_READ
Definition: IOobject.H:120