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 | Copyright (C) 2011 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 
28 #include "dictionary.H"
29 #include "volFields.H"
30 #include "volPointInterpolation.H"
32 #include "fvMesh.H"
33 #include "thresholdCellFaces.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(sampledThresholdCellFaces, 0);
41  (
42  sampledSurface,
43  sampledThresholdCellFaces,
44  word,
45  thresholdCellFaces
46  );
47 }
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 bool Foam::sampledThresholdCellFaces::updateGeometry() const
52 {
53  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
54 
55  // no update needed
56  if (fvm.time().timeIndex() == prevTimeIndex_)
57  {
58  return false;
59  }
60 
61  prevTimeIndex_ = fvm.time().timeIndex();
62 
63  // Optionally read volScalarField
64  autoPtr<volScalarField> readFieldPtr_;
65 
66  // 1. see if field in database
67  // 2. see if field can be read
68  const volScalarField* cellFldPtr = NULL;
69  if (fvm.foundObject<volScalarField>(fieldName_))
70  {
71  if (debug)
72  {
73  Info<< "sampledThresholdCellFaces::updateGeometry() : lookup "
74  << fieldName_ << endl;
75  }
76 
77  cellFldPtr = &fvm.lookupObject<volScalarField>(fieldName_);
78  }
79  else
80  {
81  // Bit of a hack. Read field and store.
82 
83  if (debug)
84  {
85  Info<< "sampledThresholdCellFaces::updateGeometry() : reading "
86  << fieldName_ << " from time " << fvm.time().timeName()
87  << endl;
88  }
89 
90  readFieldPtr_.reset
91  (
92  new volScalarField
93  (
94  IOobject
95  (
96  fieldName_,
97  fvm.time().timeName(),
98  fvm,
101  false
102  ),
103  fvm
104  )
105  );
106 
107  cellFldPtr = readFieldPtr_.operator->();
108  }
109  const volScalarField& cellFld = *cellFldPtr;
110 
111 
112  thresholdCellFaces surf
113  (
114  fvm,
115  cellFld.internalField(),
116  lowerThreshold_,
117  upperThreshold_,
118  triangulate_
119  );
120 
121  const_cast<sampledThresholdCellFaces&>
122  (
123  *this
125  meshCells_.transfer(surf.meshCells());
126 
127  // clear derived data
129 
130  if (debug)
131  {
132  Pout<< "sampledThresholdCellFaces::updateGeometry() : constructed"
133  << nl
134  << " field : " << fieldName_ << nl
135  << " lowerLimit : " << lowerThreshold_ << nl
136  << " upperLimit : " << upperThreshold_ << nl
137  << " point : " << points().size() << nl
138  << " faces : " << faces().size() << nl
139  << " cut cells : " << meshCells_.size() << endl;
140  }
141 
142  return true;
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
147 
149 (
150  const word& name,
151  const polyMesh& mesh,
152  const dictionary& dict
153 )
154 :
155  sampledSurface(name, mesh, dict),
156  fieldName_(dict.lookup("field")),
157  lowerThreshold_(dict.lookupOrDefault<scalar>("lowerLimit", -VGREAT)),
158  upperThreshold_(dict.lookupOrDefault<scalar>("upperLimit", VGREAT)),
159  zoneKey_(keyType::null),
160  triangulate_(dict.lookupOrDefault("triangulate", false)),
161  prevTimeIndex_(-1),
162  meshCells_(0)
163 {
164  if (!dict.found("lowerLimit") && !dict.found("upperLimit"))
165  {
167  (
168  "sampledThresholdCellFaces::sampledThresholdCellFaces(..)"
169  )
170  << "require at least one of 'lowerLimit' or 'upperLimit'" << endl
171  << abort(FatalError);
172  }
173 
174 // dict.readIfPresent("zone", zoneKey_);
175 //
176 // if (debug && zoneKey_.size() && mesh.cellZones().findZoneID(zoneKey_) < 0)
177 // {
178 // Info<< "cellZone " << zoneKey_
179 // << " not found - using entire mesh" << endl;
180 // }
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
185 
187 {}
188 
189 
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 
193 {
194  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
195 
196  return fvm.time().timeIndex() != prevTimeIndex_;
197 }
198 
199 
201 {
202  // already marked as expired
203  if (prevTimeIndex_ == -1)
204  {
205  return false;
206  }
207 
208  // force update
209  prevTimeIndex_ = -1;
210  return true;
211 }
212 
213 
215 {
216  return updateGeometry();
217 }
218 
219 
221 (
222  const volScalarField& vField
223 ) const
224 {
225  return sampleField(vField);
226 }
227 
228 
230 (
231  const volVectorField& vField
232 ) const
233 {
234  return sampleField(vField);
235 }
236 
237 
239 (
240  const volSphericalTensorField& vField
241 ) const
242 {
243  return sampleField(vField);
244 }
245 
246 
248 (
249  const volSymmTensorField& vField
250 ) const
251 {
252  return sampleField(vField);
253 }
254 
255 
257 (
258  const volTensorField& vField
259 ) const
260 {
261  return sampleField(vField);
262 }
263 
264 
266 (
267  const interpolation<scalar>& interpolator
268 ) const
269 {
270  return interpolateField(interpolator);
271 }
272 
273 
275 (
276  const interpolation<vector>& interpolator
277 ) const
278 {
279  return interpolateField(interpolator);
280 }
281 
284 (
285  const interpolation<sphericalTensor>& interpolator
286 ) const
287 {
288  return interpolateField(interpolator);
289 }
290 
291 
293 (
294  const interpolation<symmTensor>& interpolator
295 ) const
296 {
297  return interpolateField(interpolator);
298 }
299 
300 
302 (
303  const interpolation<tensor>& interpolator
304 ) const
305 {
306  return interpolateField(interpolator);
307 }
308 
309 
311 {
312  os << "sampledThresholdCellFaces: " << name() << " :"
313  << " field:" << fieldName_
314  << " lowerLimit:" << lowerThreshold_
315  << " upperLimit:" << upperThreshold_;
316  //<< " faces:" << faces().size() // possibly no geom yet
317  //<< " points:" << points().size();
318 }
319 
320 
321 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
const pointField & points
virtual void print(Ostream &) const
Write.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
An abstract class for surfaces with sampling.
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.
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
A class for handling words, derived from string.
Definition: word.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
messageStream Info
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
label timeIndex() const
Return current time index.
Definition: TimeState.C:73
static const keyType null
An empty keyType.
Definition: keyType.H:75
virtual bool needsUpdate() const
Does the surface need an update?
virtual void clearGeom() const
Macros for easy insertion into run-time selection tables.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
bool interpolate() const
Interpolation requested for surface.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual ~sampledThresholdCellFaces()
Destructor.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual bool update()
Update the surface as required.
void transfer(MeshedSurface< face > &)
Transfer the contents of the argument and annul the argument.
addNamedToRunTimeSelectionTable(GAMGProcAgglomeration, noneGAMGProcAgglomeration, GAMGAgglomeration, none)
sampledThresholdCellFaces(const word &name, const polyMesh &, const dictionary &)
Construct from dictionary.
A class for managing temporary objects.
Definition: PtrList.H:118
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53