sampledThresholdCellFaces.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-2021 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 
27 #include "thresholdCellFaces.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace sampledSurfaces
35 {
36  defineTypeNameAndDebug(thresholdCellFaces, 0);
37  addToRunTimeSelectionTable(sampledSurface, thresholdCellFaces, word);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 bool Foam::sampledSurfaces::thresholdCellFaces::updateGeometry() const
45 {
46  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
47 
48  // no update needed
49  if (fvm.time().timeIndex() == prevTimeIndex_)
50  {
51  return false;
52  }
53 
54  prevTimeIndex_ = fvm.time().timeIndex();
55 
56  // Optionally read volScalarField
57  autoPtr<volScalarField> readFieldPtr_;
58 
59  // 1. see if field in database
60  // 2. see if field can be read
61  const volScalarField* cellFldPtr = nullptr;
62  if (fvm.foundObject<volScalarField>(fieldName_))
63  {
64  if (debug)
65  {
66  InfoInFunction<< "Lookup " << fieldName_ << endl;
67  }
68 
69  cellFldPtr = &fvm.lookupObject<volScalarField>(fieldName_);
70  }
71  else
72  {
73  // Bit of a hack. Read field and store.
74 
75  if (debug)
76  {
78  << "Reading " << fieldName_
79  << " from time " << fvm.time().timeName()
80  << endl;
81  }
82 
83  readFieldPtr_.reset
84  (
85  new volScalarField
86  (
87  IOobject
88  (
89  fieldName_,
90  fvm.time().timeName(),
91  fvm,
94  false
95  ),
96  fvm
97  )
98  );
99 
100  cellFldPtr = readFieldPtr_.operator->();
101  }
102  const volScalarField& cellFld = *cellFldPtr;
103 
104 
106  (
107  fvm,
108  cellFld.primitiveField(),
109  lowerThreshold_,
110  upperThreshold_,
111  triangulate_
112  );
113 
114  const_cast<thresholdCellFaces&>
115  (
116  *this
118  meshCells_.transfer(surf.meshCells());
119 
120  // clear derived data
122 
123  if (debug)
124  {
125  Pout<< "thresholdCellFaces::updateGeometry() : constructed"
126  << nl
127  << " field : " << fieldName_ << nl
128  << " lowerLimit : " << lowerThreshold_ << nl
129  << " upperLimit : " << upperThreshold_ << nl
130  << " point : " << points().size() << nl
131  << " faces : " << faces().size() << nl
132  << " cut cells : " << meshCells_.size() << endl;
133  }
134 
135  return true;
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
140 
142 (
143  const word& name,
144  const polyMesh& mesh,
145  const dictionary& dict
146 )
147 :
148  sampledSurface(name, mesh, dict),
149  fieldName_(dict.lookup("field")),
150  lowerThreshold_(dict.lookupOrDefault<scalar>("lowerLimit", -vGreat)),
151  upperThreshold_(dict.lookupOrDefault<scalar>("upperLimit", vGreat)),
152  triangulate_(dict.lookupOrDefault("triangulate", false)),
153  prevTimeIndex_(-1),
154  meshCells_(0)
155 {
156  if (!dict.found("lowerLimit") && !dict.found("upperLimit"))
157  {
159  << "require at least one of 'lowerLimit' or 'upperLimit'" << endl
160  << abort(FatalError);
161  }
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
166 
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
174 {
175  return wordList(fieldName_);
176 }
177 
178 
180 {
181  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
182 
183  return fvm.time().timeIndex() != prevTimeIndex_;
184 }
185 
186 
188 {
189  // already marked as expired
190  if (prevTimeIndex_ == -1)
191  {
192  return false;
193  }
194 
195  // force update
196  prevTimeIndex_ = -1;
197  return true;
198 }
199 
200 
202 {
203  return updateGeometry();
204 }
205 
206 
209 (
210  const volScalarField& vField
211 ) const
212 {
213  return sampleField(vField);
214 }
215 
216 
219 (
220  const volVectorField& vField
221 ) const
222 {
223  return sampleField(vField);
224 }
225 
226 
229 (
230  const volSphericalTensorField& vField
231 ) const
232 {
233  return sampleField(vField);
234 }
235 
236 
239 (
240  const volSymmTensorField& vField
241 ) const
242 {
243  return sampleField(vField);
244 }
245 
246 
249 (
250  const volTensorField& vField
251 ) const
252 {
253  return sampleField(vField);
254 }
255 
256 
259 (
260  const interpolation<scalar>& interpolator
261 ) const
262 {
263  return interpolateField(interpolator);
264 }
265 
266 
269 (
270  const interpolation<vector>& interpolator
271 ) const
272 {
273  return interpolateField(interpolator);
274 }
275 
278 (
279  const interpolation<sphericalTensor>& interpolator
280 ) const
281 {
282  return interpolateField(interpolator);
283 }
284 
285 
288 (
289  const interpolation<symmTensor>& interpolator
290 ) const
291 {
292  return interpolateField(interpolator);
293 }
294 
295 
298 (
299  const interpolation<tensor>& interpolator
300 ) const
301 {
302  return interpolateField(interpolator);
303 }
304 
305 
307 {
308  os << "thresholdCellFaces: " << name() << " :"
309  << " field:" << fieldName_
310  << " lowerLimit:" << lowerThreshold_
311  << " upperLimit:" << upperThreshold_;
312  //<< " faces:" << faces().size() // possibly no geom yet
313  //<< " points:" << points().size();
314 }
315 
316 
317 // ************************************************************************* //
virtual bool update()
Update the surface as required.
Selects the mesh cell faces specified by a threshold value. Non-triangulated by default.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
thresholdCellFaces(const word &name, const polyMesh &, const dictionary &)
Construct from dictionary.
virtual wordList fields() const
Return the list of fields required.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
An abstract class for surfaces with sampling.
defineTypeNameAndDebug(distanceSurface, 0)
bool interpolate() const
Interpolation requested for surface.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void transfer(MeshedSurface< face > &)
Transfer the contents of the argument and annul the argument.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
fvMesh & mesh
virtual bool needsUpdate() const
Does the surface need an update?
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
virtual void clearGeom() const
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
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:54
static const char nl
Definition: Ostream.H:260
addToRunTimeSelectionTable(sampledSurface, distanceSurface, word)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< word > wordList
A List of words.
Definition: fileName.H:54
virtual bool expire()
Mark the surface as needing an update.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void print(Ostream &) const
Write.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
#define InfoInFunction
Report an information message using Foam::Info.