sampledThresholdCellFaces.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) 2018-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30#include "thresholdCellFaces.H"
31#include "dictionary.H"
32#include "volFields.H"
35#include "fvMesh.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
43 (
46 word,
48 );
49}
50
51
52// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53
54bool Foam::sampledThresholdCellFaces::updateGeometry() const
55{
56 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
57
58 // no update needed
59 if (fvm.time().timeIndex() == prevTimeIndex_)
60 {
61 return false;
62 }
63
64 prevTimeIndex_ = fvm.time().timeIndex();
65
66 // Use volField from database, or try to read it in
67
68 const auto* cellFldPtr = fvm.findObject<volScalarField>(fieldName_);
69
70 if (debug)
71 {
72 if (cellFldPtr)
73 {
74 InfoInFunction << "Lookup " << fieldName_ << endl;
75 }
76 else
77 {
79 << "Reading " << fieldName_
80 << " from time " << fvm.time().timeName()
81 << endl;
82 }
83 }
84
85 // For holding the volScalarField read in.
86 autoPtr<volScalarField> fieldReadPtr;
87
88 if (!cellFldPtr)
89 {
90 // Bit of a hack. Read field and store.
91
92 fieldReadPtr = autoPtr<volScalarField>::New
93 (
94 IOobject
95 (
96 fieldName_,
97 fvm.time().timeName(),
98 fvm,
101 false
102 ),
103 fvm
104 );
105 }
106
107 const volScalarField& cellFld =
108 (fieldReadPtr ? *fieldReadPtr : *cellFldPtr);
109
110 Mesh& mySurface = const_cast<sampledThresholdCellFaces&>(*this);
111
112 thresholdCellFaces surf
113 (
114 fvm,
115 cellFld.primitiveField(),
116 lowerThreshold_,
117 upperThreshold_,
118 triangulate_
119 );
120
121 mySurface.transfer(static_cast<Mesh&>(surf));
122 meshCells_.transfer(surf.meshCells());
123
124 // Clear derived data
126
127 if (debug)
128 {
129 Pout<< "sampledThresholdCellFaces::updateGeometry() : constructed"
130 << nl
131 << " field : " << fieldName_ << nl
132 << " lowerLimit : " << lowerThreshold_ << nl
133 << " upperLimit : " << upperThreshold_ << nl
134 << " point : " << points().size() << nl
135 << " faces : " << faces().size() << nl
136 << " cut cells : " << meshCells_.size() << endl;
137 }
138
139 return true;
140}
141
142
143// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
144
146(
147 const word& name,
148 const polyMesh& mesh,
149 const dictionary& dict
150)
151:
153 fieldName_(dict.get<word>("field")),
154 lowerThreshold_(dict.getOrDefault<scalar>("lowerLimit", -VGREAT)),
155 upperThreshold_(dict.getOrDefault<scalar>("upperLimit", VGREAT)),
156 triangulate_(dict.getOrDefault("triangulate", false)),
157 prevTimeIndex_(-1),
158 meshCells_()
159{
160 if (!dict.found("lowerLimit") && !dict.found("upperLimit"))
161 {
163 << "require at least one of 'lowerLimit' or 'upperLimit'" << endl
164 << abort(FatalError);
165 }
166}
167
168
169// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
170
172{
173 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
174
175 return fvm.time().timeIndex() != prevTimeIndex_;
176}
177
178
180{
181 // already marked as expired
182 if (prevTimeIndex_ == -1)
183 {
184 return false;
185 }
186
187 // force update
188 prevTimeIndex_ = -1;
189 return true;
190}
191
192
194{
195 return updateGeometry();
196}
197
198
200(
201 const interpolation<scalar>& sampler
202) const
203{
204 return sampleOnFaces(sampler);
205}
206
207
209(
210 const interpolation<vector>& sampler
211) const
212{
213 return sampleOnFaces(sampler);
214}
215
216
218(
219 const interpolation<sphericalTensor>& sampler
220) const
221{
222 return sampleOnFaces(sampler);
223}
224
225
227(
228 const interpolation<symmTensor>& sampler
229) const
230{
231 return sampleOnFaces(sampler);
232}
233
234
236(
237 const interpolation<tensor>& sampler
238) const
239{
240 return sampleOnFaces(sampler);
241}
242
243
245(
246 const interpolation<scalar>& interpolator
247) const
248{
249 return sampleOnPoints(interpolator);
250}
251
252
254(
255 const interpolation<vector>& interpolator
256) const
257{
258 return sampleOnPoints(interpolator);
259}
260
263(
264 const interpolation<sphericalTensor>& interpolator
265) const
266{
267 return sampleOnPoints(interpolator);
268}
269
270
272(
273 const interpolation<symmTensor>& interpolator
274) const
275{
276 return sampleOnPoints(interpolator);
277}
278
279
281(
282 const interpolation<tensor>& interpolator
283) const
284{
285 return sampleOnPoints(interpolator);
286}
287
288
290{
291 os << "sampledThresholdCellFaces: " << name() << " :"
292 << " field:" << fieldName_
293 << " lowerLimit:" << lowerThreshold_
294 << " upperLimit:" << upperThreshold_;
295
296 // Possibly no geom yet...
297 // if (level)
298 // {
299 // os << " faces:" << faces().size()
300 // << " points:" << points().size();
301 // }
302}
303
304
305// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
Minimal example by using system/controlDict.functions:
void transfer(List< T > &list)
Definition: List.C:447
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
label timeIndex() const noexcept
Return current time index.
Definition: TimeStateI.H:37
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
scalar print()
Print to screen.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
An abstract class for surfaces with sampling.
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
bool interpolate() const noexcept
Same as isPointData()
A sampledSurface defined by the cell faces corresponding to a threshold value.
virtual const pointField & points() const
Points of surface.
sampledThresholdCellFaces(const word &name, const polyMesh &, const dictionary &)
Construct from dictionary.
virtual const faceList & faces() const
Faces of surface.
virtual bool expire()
Mark the surface as needing an update.
virtual bool needsUpdate() const
Does the surface need an update?
virtual bool update()
Update the surface as required.
Selects the mesh cell faces specified by a threshold value. Non-triangulated by default.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict