FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

sampledThresholdCellFaces.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "sampledThresholdCellFaces.H"
00027 
00028 #include <OpenFOAM/dictionary.H>
00029 #include <finiteVolume/volFields.H>
00030 #include <finiteVolume/volPointInterpolation.H>
00031 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00032 #include <finiteVolume/fvMesh.H>
00033 #include "thresholdCellFaces.H"
00034 
00035 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00036 
00037 namespace Foam
00038 {
00039     defineTypeNameAndDebug(sampledThresholdCellFaces, 0);
00040     addNamedToRunTimeSelectionTable
00041     (
00042         sampledSurface,
00043         sampledThresholdCellFaces,
00044         word,
00045         thresholdCellFaces
00046     );
00047 }
00048 
00049 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00050 
00051 bool Foam::sampledThresholdCellFaces::updateGeometry() const
00052 {
00053     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
00054 
00055     // no update needed
00056     if (fvm.time().timeIndex() == prevTimeIndex_)
00057     {
00058         return false;
00059     }
00060 
00061     prevTimeIndex_ = fvm.time().timeIndex();
00062 
00063     // Optionally read volScalarField
00064     autoPtr<volScalarField> readFieldPtr_;
00065 
00066     // 1. see if field in database
00067     // 2. see if field can be read
00068     const volScalarField* cellFldPtr = NULL;
00069     if (fvm.foundObject<volScalarField>(fieldName_))
00070     {
00071         if (debug)
00072         {
00073             Info<< "sampledThresholdCellFaces::updateGeometry() : lookup "
00074                 << fieldName_ << endl;
00075         }
00076 
00077         cellFldPtr = &fvm.lookupObject<volScalarField>(fieldName_);
00078     }
00079     else
00080     {
00081         // Bit of a hack. Read field and store.
00082 
00083         if (debug)
00084         {
00085             Info<< "sampledThresholdCellFaces::updateGeometry() : reading "
00086                 << fieldName_ << " from time " << fvm.time().timeName()
00087                 << endl;
00088         }
00089 
00090         readFieldPtr_.reset
00091         (
00092             new volScalarField
00093             (
00094                 IOobject
00095                 (
00096                     fieldName_,
00097                     fvm.time().timeName(),
00098                     fvm,
00099                     IOobject::MUST_READ,
00100                     IOobject::NO_WRITE,
00101                     false
00102                 ),
00103                 fvm
00104             )
00105         );
00106 
00107         cellFldPtr = readFieldPtr_.operator->();
00108     }
00109     const volScalarField& cellFld = *cellFldPtr;
00110 
00111 
00112     thresholdCellFaces surf
00113     (
00114         fvm,
00115         cellFld.internalField(),
00116         lowerThreshold_,
00117         upperThreshold_,
00118         triangulate_
00119     );
00120 
00121     const_cast<sampledThresholdCellFaces&>
00122     (
00123         *this
00124     ).MeshedSurface<face>::transfer(surf);
00125     meshCells_.transfer(surf.meshCells());
00126 
00127     // clear derived data
00128     sampledSurface::clearGeom();
00129 
00130     if (debug)
00131     {
00132         Pout<< "sampledThresholdCellFaces::updateGeometry() : constructed"
00133             << nl
00134             << "    field         : " << fieldName_ << nl
00135             << "    lowerLimit    : " << lowerThreshold_ << nl
00136             << "    upperLimit    : " << upperThreshold_ << nl
00137             << "    point         : " << points().size() << nl
00138             << "    faces         : " << faces().size() << nl
00139             << "    cut cells     : " << meshCells_.size() << endl;
00140     }
00141 
00142     return true;
00143 }
00144 
00145 
00146 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00147 
00148 Foam::sampledThresholdCellFaces::sampledThresholdCellFaces
00149 (
00150     const word& name,
00151     const polyMesh& mesh,
00152     const dictionary& dict
00153 )
00154 :
00155     sampledSurface(name, mesh, dict),
00156     fieldName_(dict.lookup("field")),
00157     lowerThreshold_(dict.lookupOrDefault<scalar>("lowerLimit", -VGREAT)),
00158     upperThreshold_(dict.lookupOrDefault<scalar>("upperLimit", VGREAT)),
00159     zoneName_(word::null),
00160     triangulate_(dict.lookupOrDefault("triangulate", false)),
00161     prevTimeIndex_(-1),
00162     meshCells_(0)
00163 {
00164     if (!dict.found("lowerLimit") && !dict.found("upperLimit"))
00165     {
00166         FatalErrorIn
00167             (
00168                 "sampledThresholdCellFaces::sampledThresholdCellFaces(..)"
00169             )
00170             << "require at least one of 'lowerLimit' or 'upperLimit'" << endl
00171             << abort(FatalError);
00172     }
00173 
00174 
00175 //    dict.readIfPresent("zone", zoneName_);
00176 //
00177 //    if (debug && zoneName_.size())
00178 //    {
00179 //        if (mesh.cellZones().findZoneID(zoneName_) < 0)
00180 //        {
00181 //            Info<< "cellZone \"" << zoneName_
00182 //                << "\" not found - using entire mesh" << endl;
00183 //        }
00184 //    }
00185 }
00186 
00187 
00188 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00189 
00190 Foam::sampledThresholdCellFaces::~sampledThresholdCellFaces()
00191 {}
00192 
00193 
00194 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00195 
00196 bool Foam::sampledThresholdCellFaces::needsUpdate() const
00197 {
00198     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
00199 
00200     return fvm.time().timeIndex() != prevTimeIndex_;
00201 }
00202 
00203 
00204 bool Foam::sampledThresholdCellFaces::expire()
00205 {
00206     // already marked as expired
00207     if (prevTimeIndex_ == -1)
00208     {
00209         return false;
00210     }
00211 
00212     // force update
00213     prevTimeIndex_ = -1;
00214     return true;
00215 }
00216 
00217 
00218 bool Foam::sampledThresholdCellFaces::update()
00219 {
00220     return updateGeometry();
00221 }
00222 
00223 
00224 Foam::tmp<Foam::scalarField>
00225 Foam::sampledThresholdCellFaces::sample
00226 (
00227     const volScalarField& vField
00228 ) const
00229 {
00230     return sampleField(vField);
00231 }
00232 
00233 
00234 Foam::tmp<Foam::vectorField>
00235 Foam::sampledThresholdCellFaces::sample
00236 (
00237     const volVectorField& vField
00238 ) const
00239 {
00240     return sampleField(vField);
00241 }
00242 
00243 
00244 Foam::tmp<Foam::sphericalTensorField>
00245 Foam::sampledThresholdCellFaces::sample
00246 (
00247     const volSphericalTensorField& vField
00248 ) const
00249 {
00250     return sampleField(vField);
00251 }
00252 
00253 
00254 Foam::tmp<Foam::symmTensorField>
00255 Foam::sampledThresholdCellFaces::sample
00256 (
00257     const volSymmTensorField& vField
00258 ) const
00259 {
00260     return sampleField(vField);
00261 }
00262 
00263 
00264 Foam::tmp<Foam::tensorField>
00265 Foam::sampledThresholdCellFaces::sample
00266 (
00267     const volTensorField& vField
00268 ) const
00269 {
00270     return sampleField(vField);
00271 }
00272 
00273 
00274 Foam::tmp<Foam::scalarField>
00275 Foam::sampledThresholdCellFaces::interpolate
00276 (
00277     const interpolation<scalar>& interpolator
00278 ) const
00279 {
00280     return interpolateField(interpolator);
00281 }
00282 
00283 
00284 Foam::tmp<Foam::vectorField>
00285 Foam::sampledThresholdCellFaces::interpolate
00286 (
00287     const interpolation<vector>& interpolator
00288 ) const
00289 {
00290     return interpolateField(interpolator);
00291 }
00292 
00293 Foam::tmp<Foam::sphericalTensorField>
00294 Foam::sampledThresholdCellFaces::interpolate
00295 (
00296     const interpolation<sphericalTensor>& interpolator
00297 ) const
00298 {
00299     return interpolateField(interpolator);
00300 }
00301 
00302 
00303 Foam::tmp<Foam::symmTensorField>
00304 Foam::sampledThresholdCellFaces::interpolate
00305 (
00306     const interpolation<symmTensor>& interpolator
00307 ) const
00308 {
00309     return interpolateField(interpolator);
00310 }
00311 
00312 
00313 Foam::tmp<Foam::tensorField>
00314 Foam::sampledThresholdCellFaces::interpolate
00315 (
00316     const interpolation<tensor>& interpolator
00317 ) const
00318 {
00319     return interpolateField(interpolator);
00320 }
00321 
00322 
00323 void Foam::sampledThresholdCellFaces::print(Ostream& os) const
00324 {
00325     os  << "sampledThresholdCellFaces: " << name() << " :"
00326         << "  field:" << fieldName_
00327         << "  lowerLimit:" << lowerThreshold_
00328         << "  upperLimit:" << upperThreshold_;
00329         //<< "  faces:" << faces().size()   // possibly no geom yet
00330         //<< "  points:" << points().size();
00331 }
00332 
00333 
00334 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines