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-2020 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  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
176 
177  return fvm.time().timeIndex() != prevTimeIndex_;
178 }
179 
180 
182 {
183  // already marked as expired
184  if (prevTimeIndex_ == -1)
185  {
186  return false;
187  }
188 
189  // force update
190  prevTimeIndex_ = -1;
191  return true;
192 }
193 
194 
196 {
197  return updateGeometry();
198 }
199 
200 
203 (
204  const volScalarField& vField
205 ) const
206 {
207  return sampleField(vField);
208 }
209 
210 
213 (
214  const volVectorField& vField
215 ) const
216 {
217  return sampleField(vField);
218 }
219 
220 
223 (
224  const volSphericalTensorField& vField
225 ) const
226 {
227  return sampleField(vField);
228 }
229 
230 
233 (
234  const volSymmTensorField& vField
235 ) const
236 {
237  return sampleField(vField);
238 }
239 
240 
243 (
244  const volTensorField& vField
245 ) const
246 {
247  return sampleField(vField);
248 }
249 
250 
253 (
254  const interpolation<scalar>& interpolator
255 ) const
256 {
257  return interpolateField(interpolator);
258 }
259 
260 
263 (
264  const interpolation<vector>& interpolator
265 ) const
266 {
267  return interpolateField(interpolator);
268 }
269 
272 (
273  const interpolation<sphericalTensor>& interpolator
274 ) const
275 {
276  return interpolateField(interpolator);
277 }
278 
279 
282 (
283  const interpolation<symmTensor>& interpolator
284 ) const
285 {
286  return interpolateField(interpolator);
287 }
288 
289 
292 (
293  const interpolation<tensor>& interpolator
294 ) const
295 {
296  return interpolateField(interpolator);
297 }
298 
299 
301 {
302  os << "thresholdCellFaces: " << name() << " :"
303  << " field:" << fieldName_
304  << " lowerLimit:" << lowerThreshold_
305  << " upperLimit:" << upperThreshold_;
306  //<< " faces:" << faces().size() // possibly no geom yet
307  //<< " points:" << points().size();
308 }
309 
310 
311 // ************************************************************************* //
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:643
thresholdCellFaces(const word &name, const polyMesh &, const dictionary &)
Construct from dictionary.
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:323
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:239
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:57
dynamicFvMesh & mesh
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
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
addToRunTimeSelectionTable(sampledSurface, distanceSurface, word)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual bool expire()
Mark the surface as needing an update.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
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:78
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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:844
#define InfoInFunction
Report an information message using Foam::Info.