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-2015 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 
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;
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(),
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 :
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  << "require at least one of 'lowerLimit' or 'upperLimit'" << endl
168  << abort(FatalError);
169  }
170 
171 // dict.readIfPresent("zone", zoneKey_);
172 //
173 // if (debug && zoneKey_.size() && mesh.cellZones().findZoneID(zoneKey_) < 0)
174 // {
175 // Info<< "cellZone " << zoneKey_
176 // << " not found - using entire mesh" << endl;
177 // }
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
182 
184 {}
185 
186 
187 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
188 
190 {
191  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
192 
193  return fvm.time().timeIndex() != prevTimeIndex_;
194 }
195 
196 
198 {
199  // already marked as expired
200  if (prevTimeIndex_ == -1)
201  {
202  return false;
203  }
204 
205  // force update
206  prevTimeIndex_ = -1;
207  return true;
208 }
209 
210 
212 {
213  return updateGeometry();
214 }
215 
216 
218 (
219  const volScalarField& vField
220 ) const
221 {
222  return sampleField(vField);
223 }
224 
225 
227 (
228  const volVectorField& vField
229 ) const
230 {
231  return sampleField(vField);
232 }
233 
234 
236 (
237  const volSphericalTensorField& vField
238 ) const
239 {
240  return sampleField(vField);
241 }
242 
243 
245 (
246  const volSymmTensorField& vField
247 ) const
248 {
249  return sampleField(vField);
250 }
251 
252 
254 (
255  const volTensorField& vField
256 ) const
257 {
258  return sampleField(vField);
259 }
260 
261 
263 (
264  const interpolation<scalar>& interpolator
265 ) const
266 {
267  return interpolateField(interpolator);
268 }
269 
270 
272 (
273  const interpolation<vector>& interpolator
274 ) const
275 {
276  return interpolateField(interpolator);
277 }
278 
281 (
282  const interpolation<sphericalTensor>& interpolator
283 ) const
284 {
285  return interpolateField(interpolator);
286 }
287 
288 
290 (
291  const interpolation<symmTensor>& interpolator
292 ) const
293 {
294  return interpolateField(interpolator);
295 }
296 
297 
299 (
300  const interpolation<tensor>& interpolator
301 ) const
302 {
303  return interpolateField(interpolator);
304 }
305 
306 
308 {
309  os << "sampledThresholdCellFaces: " << name() << " :"
310  << " field:" << fieldName_
311  << " lowerLimit:" << lowerThreshold_
312  << " upperLimit:" << upperThreshold_;
313  //<< " faces:" << faces().size() // possibly no geom yet
314  //<< " points:" << points().size();
315 }
316 
317 
318 // ************************************************************************* //
volFields.H
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(fvPatch, cyclicAMIFvPatch, polyPatch, cyclicPeriodicAMI)
Foam::sampledThresholdCellFaces::lowerThreshold_
const scalar lowerThreshold_
Threshold value.
Definition: sampledThresholdCellFaces.H:64
Foam::sampledThresholdCellFaces::~sampledThresholdCellFaces
virtual ~sampledThresholdCellFaces()
Destructor.
Definition: sampledThresholdCellFaces.C:183
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::sampledThresholdCellFaces::prevTimeIndex_
label prevTimeIndex_
Time at last call, also track it surface needs an update.
Definition: sampledThresholdCellFaces.H:78
Foam::List::transfer
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Foam::sampledThresholdCellFaces::fieldName_
const word fieldName_
Field to get isoSurface of.
Definition: sampledThresholdCellFaces.H:61
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::sampledThresholdCellFaces::sampledThresholdCellFaces
sampledThresholdCellFaces(const word &name, const polyMesh &, const dictionary &)
Construct from dictionary.
Definition: sampledThresholdCellFaces.C:149
Foam::sampledThresholdCellFaces::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledThresholdCellFaces.C:189
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::sampledThresholdCellFaces::update
virtual bool update()
Update the surface as required.
Definition: sampledThresholdCellFaces.C:211
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::sampledThresholdCellFaces::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledThresholdCellFaces.C:197
Foam::sampledThresholdCellFaces::sample
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledThresholdCellFaces.C:218
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::sampledThresholdCellFaces::meshCells_
labelList meshCells_
For every face the original cell in mesh.
Definition: sampledThresholdCellFaces.H:81
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::thresholdCellFaces::meshCells
labelList & meshCells()
For every face original cell in mesh.
Definition: thresholdCellFaces.H:96
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:77
Foam::keyType::null
static const keyType null
An empty keyType.
Definition: keyType.H:75
thresholdCellFaces.H
Foam::interpolation< scalar >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::thresholdCellFaces
Selects the mesh cell faces specified by a threshold value. Non-triangulated by default.
Definition: thresholdCellFaces.H:51
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::sampledThresholdCellFaces::upperThreshold_
const scalar upperThreshold_
Threshold value.
Definition: sampledThresholdCellFaces.H:67
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
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:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
volPointInterpolation.H
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::sampledThresholdCellFaces
A sampledSurface defined by the cell faces corresponding to a threshold value.
Definition: sampledThresholdCellFaces.H:50
Foam::sampledThresholdCellFaces::updateGeometry
bool updateGeometry() const
Create surface (if time has changed)
Definition: sampledThresholdCellFaces.C:51
dictionary.H
Foam::sampledSurface::clearGeom
virtual void clearGeom() const
Definition: sampledSurface.C:42
Foam::sampledThresholdCellFaces::faces
virtual const faceList & faces() const
Faces of surface.
Definition: sampledThresholdCellFaces.H:145
sampledThresholdCellFaces.H
Foam::autoPtr::reset
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::sampledThresholdCellFaces::triangulate_
bool triangulate_
Triangulated faces or keep faces as is.
Definition: sampledThresholdCellFaces.H:73
Foam::sampledSurface::mesh
const polyMesh & mesh() const
Access to the underlying mesh.
Definition: sampledSurface.H:244
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::TimeState::timeIndex
label timeIndex() const
Return current time index.
Definition: TimeState.C:73
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::MeshedSurface
A surface geometry mesh with zone information, not to be confused with the similarly named surfaceMes...
Definition: MeshedSurface.H:72
Foam::sampledSurface::interpolate
bool interpolate() const
Interpolation requested for surface.
Definition: sampledSurface.H:256
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::sampledThresholdCellFaces::print
virtual void print(Ostream &) const
Write.
Definition: sampledThresholdCellFaces.C:307
Foam::sampledThresholdCellFaces::points
virtual const pointField & points() const
Points of surface.
Definition: sampledThresholdCellFaces.H:139
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47