initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / sampling / sampledSurface / thresholdCellFaces / sampledThresholdCellFaces.C
blob1dd43b4be3b899c839339e2d5c2d1970eeb8a85d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "sampledThresholdCellFaces.H"
29 #include "dictionary.H"
30 #include "volFields.H"
31 #include "volPointInterpolation.H"
32 #include "addToRunTimeSelectionTable.H"
33 #include "fvMesh.H"
34 #include "thresholdCellFaces.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 namespace Foam
40     defineTypeNameAndDebug(sampledThresholdCellFaces, 0);
41     addNamedToRunTimeSelectionTable
42     (
43         sampledSurface,
44         sampledThresholdCellFaces,
45         word,
46         thresholdCellFaces
47     );
50 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
52 bool Foam::sampledThresholdCellFaces::updateGeometry() const
54     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
56     // no update needed
57     if (fvm.time().timeIndex() == prevTimeIndex_)
58     {
59         return false;
60     }
62     prevTimeIndex_ = fvm.time().timeIndex();
64     // Optionally read volScalarField
65     autoPtr<volScalarField> readFieldPtr_;
67     // 1. see if field in database
68     // 2. see if field can be read
69     const volScalarField* cellFldPtr = NULL;
70     if (fvm.foundObject<volScalarField>(fieldName_))
71     {
72         if (debug)
73         {
74             Info<< "sampledThresholdCellFaces::updateGeometry() : lookup "
75                 << fieldName_ << endl;
76         }
78         cellFldPtr = &fvm.lookupObject<volScalarField>(fieldName_);
79     }
80     else
81     {
82         // Bit of a hack. Read field and store.
84         if (debug)
85         {
86             Info<< "sampledThresholdCellFaces::updateGeometry() : reading "
87                 << fieldName_ << " from time " << fvm.time().timeName()
88                 << endl;
89         }
91         readFieldPtr_.reset
92         (
93             new volScalarField
94             (
95                 IOobject
96                 (
97                     fieldName_,
98                     fvm.time().timeName(),
99                     fvm,
100                     IOobject::MUST_READ,
101                     IOobject::NO_WRITE,
102                     false
103                 ),
104                 fvm
105             )
106         );
108         cellFldPtr = readFieldPtr_.operator->();
109     }
110     const volScalarField& cellFld = *cellFldPtr;
113     thresholdCellFaces surf
114     (
115         fvm,
116         cellFld.internalField(),
117         lowerThreshold_,
118         upperThreshold_,
119         triangulate_
120     );
122     const_cast<sampledThresholdCellFaces&>
123     (
124         *this
125     ).MeshedSurface<face>::transfer(surf);
126     meshCells_.transfer(surf.meshCells());
129     if (debug)
130     {
131         Pout<< "sampledThresholdCellFaces::updateGeometry() : constructed"
132             << nl
133             << "    field         : " << fieldName_ << nl
134             << "    lowerLimit    : " << lowerThreshold_ << nl
135             << "    upperLimit    : " << upperThreshold_ << nl
136             << "    point         : " << points().size() << nl
137             << "    faces         : " << faces().size() << nl
138             << "    cut cells     : " << meshCells_.size() << endl;
139     }
141     return true;
145 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
147 Foam::sampledThresholdCellFaces::sampledThresholdCellFaces
149     const word& name,
150     const polyMesh& mesh,
151     const dictionary& dict
154     sampledSurface(name, mesh, dict),
155     fieldName_(dict.lookup("field")),
156     lowerThreshold_(dict.lookupOrDefault<scalar>("lowerLimit", -VGREAT)),
157     upperThreshold_(dict.lookupOrDefault<scalar>("upperLimit", VGREAT)),
158     zoneName_(word::null),
159     triangulate_(dict.lookupOrDefault("triangulate", false)),
160     prevTimeIndex_(-1),
161     meshCells_(0)
163     if (!dict.found("lowerLimit") && !dict.found("upperLimit"))
164     {
165         FatalErrorIn
166             (
167                 "sampledThresholdCellFaces::sampledThresholdCellFaces(..)"
168             )
169             << "require at least one of 'lowerLimit' or 'upperLimit'" << endl
170             << abort(FatalError);
171     }
174 //    dict.readIfPresent("zone", zoneName_);
176 //    if (debug && zoneName_.size())
177 //    {
178 //        if (mesh.cellZones().findZoneID(zoneName_) < 0)
179 //        {
180 //            Info<< "cellZone \"" << zoneName_
181 //                << "\" not found - using entire mesh" << endl;
182 //        }
183 //    }
187 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
189 Foam::sampledThresholdCellFaces::~sampledThresholdCellFaces()
193 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
195 bool Foam::sampledThresholdCellFaces::needsUpdate() const
197     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
199     return fvm.time().timeIndex() != prevTimeIndex_;
203 bool Foam::sampledThresholdCellFaces::expire()
205     // already marked as expired
206     if (prevTimeIndex_ == -1)
207     {
208         return false;
209     }
211     // force update
212     prevTimeIndex_ = -1;
213     return true;
217 bool Foam::sampledThresholdCellFaces::update()
219     return updateGeometry();
223 Foam::tmp<Foam::scalarField>
224 Foam::sampledThresholdCellFaces::sample
226     const volScalarField& vField
227 ) const
229     return sampleField(vField);
233 Foam::tmp<Foam::vectorField>
234 Foam::sampledThresholdCellFaces::sample
236     const volVectorField& vField
237 ) const
239     return sampleField(vField);
243 Foam::tmp<Foam::sphericalTensorField>
244 Foam::sampledThresholdCellFaces::sample
246     const volSphericalTensorField& vField
247 ) const
249     return sampleField(vField);
253 Foam::tmp<Foam::symmTensorField>
254 Foam::sampledThresholdCellFaces::sample
256     const volSymmTensorField& vField
257 ) const
259     return sampleField(vField);
263 Foam::tmp<Foam::tensorField>
264 Foam::sampledThresholdCellFaces::sample
266     const volTensorField& vField
267 ) const
269     return sampleField(vField);
273 Foam::tmp<Foam::scalarField>
274 Foam::sampledThresholdCellFaces::interpolate
276     const interpolation<scalar>& interpolator
277 ) const
279     return interpolateField(interpolator);
283 Foam::tmp<Foam::vectorField>
284 Foam::sampledThresholdCellFaces::interpolate
286     const interpolation<vector>& interpolator
287 ) const
289     return interpolateField(interpolator);
292 Foam::tmp<Foam::sphericalTensorField>
293 Foam::sampledThresholdCellFaces::interpolate
295     const interpolation<sphericalTensor>& interpolator
296 ) const
298     return interpolateField(interpolator);
302 Foam::tmp<Foam::symmTensorField>
303 Foam::sampledThresholdCellFaces::interpolate
305     const interpolation<symmTensor>& interpolator
306 ) const
308     return interpolateField(interpolator);
312 Foam::tmp<Foam::tensorField>
313 Foam::sampledThresholdCellFaces::interpolate
315     const interpolation<tensor>& interpolator
316 ) const
318     return interpolateField(interpolator);
322 void Foam::sampledThresholdCellFaces::print(Ostream& os) const
324     os  << "sampledThresholdCellFaces: " << name() << " :"
325         << "  field:" << fieldName_
326         << "  lowerLimit:" << lowerThreshold_
327         << "  upperLimit:" << upperThreshold_;
328         //<< "  faces:" << faces().size()   // possibly no geom yet
329         //<< "  points:" << points().size();
333 // ************************************************************************* //