1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "sampledThresholdCellFaces.H"
28 #include "dictionary.H"
29 #include "volFields.H"
30 #include "volPointInterpolation.H"
31 #include "addToRunTimeSelectionTable.H"
33 #include "thresholdCellFaces.H"
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 defineTypeNameAndDebug(sampledThresholdCellFaces, 0);
40 addNamedToRunTimeSelectionTable
43 sampledThresholdCellFaces,
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 bool Foam::sampledThresholdCellFaces::updateGeometry() const
53 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
56 if (fvm.time().timeIndex() == prevTimeIndex_)
61 prevTimeIndex_ = fvm.time().timeIndex();
63 // Optionally read volScalarField
64 autoPtr<volScalarField> readFieldPtr_;
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_))
73 Info<< "sampledThresholdCellFaces::updateGeometry() : lookup "
74 << fieldName_ << endl;
77 cellFldPtr = &fvm.lookupObject<volScalarField>(fieldName_);
81 // Bit of a hack. Read field and store.
85 Info<< "sampledThresholdCellFaces::updateGeometry() : reading "
86 << fieldName_ << " from time " << fvm.time().timeName()
97 fvm.time().timeName(),
107 cellFldPtr = readFieldPtr_.operator->();
109 const volScalarField& cellFld = *cellFldPtr;
112 thresholdCellFaces surf
115 cellFld.internalField(),
121 const_cast<sampledThresholdCellFaces&>
124 ).MeshedSurface<face>::transfer(surf);
125 meshCells_.transfer(surf.meshCells());
127 // clear derived data
128 sampledSurface::clearGeom();
132 Pout<< "sampledThresholdCellFaces::updateGeometry() : constructed"
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;
146 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
148 Foam::sampledThresholdCellFaces::sampledThresholdCellFaces
151 const polyMesh& mesh,
152 const dictionary& dict
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 zoneName_(word::null),
160 triangulate_(dict.lookupOrDefault("triangulate", false)),
164 if (!dict.found("lowerLimit") && !dict.found("upperLimit"))
168 "sampledThresholdCellFaces::sampledThresholdCellFaces(..)"
170 << "require at least one of 'lowerLimit' or 'upperLimit'" << endl
171 << abort(FatalError);
175 // dict.readIfPresent("zone", zoneName_);
177 // if (debug && zoneName_.size())
179 // if (mesh.cellZones().findZoneID(zoneName_) < 0)
181 // Info<< "cellZone \"" << zoneName_
182 // << "\" not found - using entire mesh" << endl;
188 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
190 Foam::sampledThresholdCellFaces::~sampledThresholdCellFaces()
194 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 bool Foam::sampledThresholdCellFaces::needsUpdate() const
198 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
200 return fvm.time().timeIndex() != prevTimeIndex_;
204 bool Foam::sampledThresholdCellFaces::expire()
206 // already marked as expired
207 if (prevTimeIndex_ == -1)
218 bool Foam::sampledThresholdCellFaces::update()
220 return updateGeometry();
224 Foam::tmp<Foam::scalarField>
225 Foam::sampledThresholdCellFaces::sample
227 const volScalarField& vField
230 return sampleField(vField);
234 Foam::tmp<Foam::vectorField>
235 Foam::sampledThresholdCellFaces::sample
237 const volVectorField& vField
240 return sampleField(vField);
244 Foam::tmp<Foam::sphericalTensorField>
245 Foam::sampledThresholdCellFaces::sample
247 const volSphericalTensorField& vField
250 return sampleField(vField);
254 Foam::tmp<Foam::symmTensorField>
255 Foam::sampledThresholdCellFaces::sample
257 const volSymmTensorField& vField
260 return sampleField(vField);
264 Foam::tmp<Foam::tensorField>
265 Foam::sampledThresholdCellFaces::sample
267 const volTensorField& vField
270 return sampleField(vField);
274 Foam::tmp<Foam::scalarField>
275 Foam::sampledThresholdCellFaces::interpolate
277 const interpolation<scalar>& interpolator
280 return interpolateField(interpolator);
284 Foam::tmp<Foam::vectorField>
285 Foam::sampledThresholdCellFaces::interpolate
287 const interpolation<vector>& interpolator
290 return interpolateField(interpolator);
293 Foam::tmp<Foam::sphericalTensorField>
294 Foam::sampledThresholdCellFaces::interpolate
296 const interpolation<sphericalTensor>& interpolator
299 return interpolateField(interpolator);
303 Foam::tmp<Foam::symmTensorField>
304 Foam::sampledThresholdCellFaces::interpolate
306 const interpolation<symmTensor>& interpolator
309 return interpolateField(interpolator);
313 Foam::tmp<Foam::tensorField>
314 Foam::sampledThresholdCellFaces::interpolate
316 const interpolation<tensor>& interpolator
319 return interpolateField(interpolator);
323 void Foam::sampledThresholdCellFaces::print(Ostream& os) const
325 os << "sampledThresholdCellFaces: " << name() << " :"
326 << " field:" << fieldName_
327 << " lowerLimit:" << lowerThreshold_
328 << " upperLimit:" << upperThreshold_;
329 //<< " faces:" << faces().size() // possibly no geom yet
330 //<< " points:" << points().size();
334 // ************************************************************************* //