ENH: patchCloudSet: new sampledSet
[OpenFOAM-1.7.x.git] / src / sampling / sampledSurface / thresholdCellFaces / sampledThresholdCellFaces.C
blobc6f2642b6fe0d33f942195045e8129cc36d798d4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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"
32 #include "fvMesh.H"
33 #include "thresholdCellFaces.H"
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 namespace Foam
39     defineTypeNameAndDebug(sampledThresholdCellFaces, 0);
40     addNamedToRunTimeSelectionTable
41     (
42         sampledSurface,
43         sampledThresholdCellFaces,
44         word,
45         thresholdCellFaces
46     );
49 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
51 bool Foam::sampledThresholdCellFaces::updateGeometry() const
53     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
55     // no update needed
56     if (fvm.time().timeIndex() == prevTimeIndex_)
57     {
58         return false;
59     }
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_))
70     {
71         if (debug)
72         {
73             Info<< "sampledThresholdCellFaces::updateGeometry() : lookup "
74                 << fieldName_ << endl;
75         }
77         cellFldPtr = &fvm.lookupObject<volScalarField>(fieldName_);
78     }
79     else
80     {
81         // Bit of a hack. Read field and store.
83         if (debug)
84         {
85             Info<< "sampledThresholdCellFaces::updateGeometry() : reading "
86                 << fieldName_ << " from time " << fvm.time().timeName()
87                 << endl;
88         }
90         readFieldPtr_.reset
91         (
92             new volScalarField
93             (
94                 IOobject
95                 (
96                     fieldName_,
97                     fvm.time().timeName(),
98                     fvm,
99                     IOobject::MUST_READ,
100                     IOobject::NO_WRITE,
101                     false
102                 ),
103                 fvm
104             )
105         );
107         cellFldPtr = readFieldPtr_.operator->();
108     }
109     const volScalarField& cellFld = *cellFldPtr;
112     thresholdCellFaces surf
113     (
114         fvm,
115         cellFld.internalField(),
116         lowerThreshold_,
117         upperThreshold_,
118         triangulate_
119     );
121     const_cast<sampledThresholdCellFaces&>
122     (
123         *this
124     ).MeshedSurface<face>::transfer(surf);
125     meshCells_.transfer(surf.meshCells());
127     // clear derived data
128     sampledSurface::clearGeom();
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     }
142     return true;
146 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
148 Foam::sampledThresholdCellFaces::sampledThresholdCellFaces
150     const word& name,
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)),
161     prevTimeIndex_(-1),
162     meshCells_(0)
164     if (!dict.found("lowerLimit") && !dict.found("upperLimit"))
165     {
166         FatalErrorIn
167             (
168                 "sampledThresholdCellFaces::sampledThresholdCellFaces(..)"
169             )
170             << "require at least one of 'lowerLimit' or 'upperLimit'" << endl
171             << abort(FatalError);
172     }
175 //    dict.readIfPresent("zone", zoneName_);
177 //    if (debug && zoneName_.size())
178 //    {
179 //        if (mesh.cellZones().findZoneID(zoneName_) < 0)
180 //        {
181 //            Info<< "cellZone \"" << zoneName_
182 //                << "\" not found - using entire mesh" << endl;
183 //        }
184 //    }
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)
208     {
209         return false;
210     }
212     // force update
213     prevTimeIndex_ = -1;
214     return true;
218 bool Foam::sampledThresholdCellFaces::update()
220     return updateGeometry();
224 Foam::tmp<Foam::scalarField>
225 Foam::sampledThresholdCellFaces::sample
227     const volScalarField& vField
228 ) const
230     return sampleField(vField);
234 Foam::tmp<Foam::vectorField>
235 Foam::sampledThresholdCellFaces::sample
237     const volVectorField& vField
238 ) const
240     return sampleField(vField);
244 Foam::tmp<Foam::sphericalTensorField>
245 Foam::sampledThresholdCellFaces::sample
247     const volSphericalTensorField& vField
248 ) const
250     return sampleField(vField);
254 Foam::tmp<Foam::symmTensorField>
255 Foam::sampledThresholdCellFaces::sample
257     const volSymmTensorField& vField
258 ) const
260     return sampleField(vField);
264 Foam::tmp<Foam::tensorField>
265 Foam::sampledThresholdCellFaces::sample
267     const volTensorField& vField
268 ) const
270     return sampleField(vField);
274 Foam::tmp<Foam::scalarField>
275 Foam::sampledThresholdCellFaces::interpolate
277     const interpolation<scalar>& interpolator
278 ) const
280     return interpolateField(interpolator);
284 Foam::tmp<Foam::vectorField>
285 Foam::sampledThresholdCellFaces::interpolate
287     const interpolation<vector>& interpolator
288 ) const
290     return interpolateField(interpolator);
293 Foam::tmp<Foam::sphericalTensorField>
294 Foam::sampledThresholdCellFaces::interpolate
296     const interpolation<sphericalTensor>& interpolator
297 ) const
299     return interpolateField(interpolator);
303 Foam::tmp<Foam::symmTensorField>
304 Foam::sampledThresholdCellFaces::interpolate
306     const interpolation<symmTensor>& interpolator
307 ) const
309     return interpolateField(interpolator);
313 Foam::tmp<Foam::tensorField>
314 Foam::sampledThresholdCellFaces::interpolate
316     const interpolation<tensor>& interpolator
317 ) const
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 // ************************************************************************* //