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

sampledIsoSurfaceCell.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 "sampledIsoSurfaceCell.H"
00027 #include <OpenFOAM/dictionary.H>
00028 #include <finiteVolume/volFields.H>
00029 #include <finiteVolume/volPointInterpolation.H>
00030 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00031 #include <finiteVolume/fvMesh.H>
00032 #include "isoSurfaceCell.H"
00033 
00034 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00035 
00036 namespace Foam
00037 {
00038     defineTypeNameAndDebug(sampledIsoSurfaceCell, 0);
00039     addNamedToRunTimeSelectionTable(sampledSurface, sampledIsoSurfaceCell, word, isoSurfaceCell);
00040 }
00041 
00042 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00043 
00044 bool Foam::sampledIsoSurfaceCell::updateGeometry() const
00045 {
00046     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
00047 
00048     // no update needed
00049     if (fvm.time().timeIndex() == prevTimeIndex_)
00050     {
00051         return false;
00052     }
00053 
00054     prevTimeIndex_ = fvm.time().timeIndex();
00055 
00056     // Clear any stored topo
00057     facesPtr_.clear();
00058 
00059     // Clear derived data
00060     sampledSurface::clearGeom();
00061 
00062     // Optionally read volScalarField
00063     autoPtr<volScalarField> readFieldPtr_;
00064 
00065     // 1. see if field in database
00066     // 2. see if field can be read
00067     const volScalarField* cellFldPtr = NULL;
00068     if (fvm.foundObject<volScalarField>(isoField_))
00069     {
00070         if (debug)
00071         {
00072             Info<< "sampledIsoSurfaceCell::updateGeometry() : lookup "
00073                 << isoField_ << endl;
00074         }
00075 
00076         cellFldPtr = &fvm.lookupObject<volScalarField>(isoField_);
00077     }
00078     else
00079     {
00080         // Bit of a hack. Read field and store.
00081 
00082         if (debug)
00083         {
00084             Info<< "sampledIsoSurfaceCell::updateGeometry() : reading "
00085                 << isoField_ << " from time " <<fvm.time().timeName()
00086                 << endl;
00087         }
00088 
00089         readFieldPtr_.reset
00090         (
00091             new volScalarField
00092             (
00093                 IOobject
00094                 (
00095                     isoField_,
00096                     fvm.time().timeName(),
00097                     fvm,
00098                     IOobject::MUST_READ,
00099                     IOobject::NO_WRITE,
00100                     false
00101                 ),
00102                 fvm
00103             )
00104         );
00105 
00106         cellFldPtr = readFieldPtr_.operator->();
00107     }
00108     const volScalarField& cellFld = *cellFldPtr;
00109 
00110     tmp<pointScalarField> pointFld
00111     (
00112         volPointInterpolation::New(fvm).interpolate(cellFld)
00113     );
00114 
00115     if (average_)
00116     {
00117         //- From point field and interpolated cell.
00118         scalarField cellAvg(fvm.nCells(), scalar(0.0));
00119         labelField nPointCells(fvm.nCells(), 0);
00120         {
00121             for (label pointI = 0; pointI < fvm.nPoints(); pointI++)
00122             {
00123                 const labelList& pCells = fvm.pointCells(pointI);
00124 
00125                 forAll(pCells, i)
00126                 {
00127                     label cellI = pCells[i];
00128 
00129                     cellAvg[cellI] += pointFld().internalField()[pointI];
00130                     nPointCells[cellI]++;
00131                 }
00132             }
00133         }
00134         forAll(cellAvg, cellI)
00135         {
00136             cellAvg[cellI] /= nPointCells[cellI];
00137         }
00138 
00139         const isoSurfaceCell iso
00140         (
00141             fvm,
00142             cellAvg,
00143             pointFld().internalField(),
00144             isoVal_,
00145             regularise_
00146         );
00147 
00148         const_cast<sampledIsoSurfaceCell&>
00149         (
00150             *this
00151         ).triSurface::operator=(iso);
00152         meshCells_ = iso.meshCells();
00153     }
00154     else
00155     {
00156         //- Direct from cell field and point field. Gives bad continuity.
00157         const isoSurfaceCell iso
00158         (
00159             fvm,
00160             cellFld.internalField(),
00161             pointFld().internalField(),
00162             isoVal_,
00163             regularise_
00164         );
00165 
00166         const_cast<sampledIsoSurfaceCell&>
00167         (
00168             *this
00169         ).triSurface::operator=(iso);
00170         meshCells_ = iso.meshCells();
00171     }
00172 
00173 
00174     if (debug)
00175     {
00176         Pout<< "sampledIsoSurfaceCell::updateGeometry() : constructed iso:"
00177             << nl
00178             << "    regularise     : " << regularise_ << nl
00179             << "    average        : " << average_ << nl
00180             << "    isoField       : " << isoField_ << nl
00181             << "    isoValue       : " << isoVal_ << nl
00182             << "    points         : " << points().size() << nl
00183             << "    tris           : " << triSurface::size() << nl
00184             << "    cut cells      : " << meshCells_.size() << endl;
00185     }
00186 
00187     return true;
00188 }
00189 
00190 
00191 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00192 
00193 Foam::sampledIsoSurfaceCell::sampledIsoSurfaceCell
00194 (
00195     const word& name,
00196     const polyMesh& mesh,
00197     const dictionary& dict
00198 )
00199 :
00200     sampledSurface(name, mesh, dict),
00201     isoField_(dict.lookup("isoField")),
00202     isoVal_(readScalar(dict.lookup("isoValue"))),
00203     regularise_(dict.lookupOrDefault("regularise", true)),
00204     average_(dict.lookupOrDefault("average", true)),
00205     zoneName_(word::null),
00206     facesPtr_(NULL),
00207     prevTimeIndex_(-1),
00208     meshCells_(0)
00209 {
00210 //    dict.readIfPresent("zone", zoneName_);
00211 //
00212 //    if (debug && zoneName_.size())
00213 //    {
00214 //        if (mesh.cellZones().findZoneID(zoneName_) < 0)
00215 //        {
00216 //            Info<< "cellZone \"" << zoneName_
00217 //                << "\" not found - using entire mesh" << endl;
00218 //        }
00219 //    }
00220 }
00221 
00222 
00223 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00224 
00225 Foam::sampledIsoSurfaceCell::~sampledIsoSurfaceCell()
00226 {}
00227 
00228 
00229 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00230 
00231 bool Foam::sampledIsoSurfaceCell::needsUpdate() const
00232 {
00233     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
00234 
00235     return fvm.time().timeIndex() != prevTimeIndex_;
00236 }
00237 
00238 
00239 bool Foam::sampledIsoSurfaceCell::expire()
00240 {
00241     facesPtr_.clear();
00242 
00243     // Clear derived data
00244     sampledSurface::clearGeom();
00245 
00246     // already marked as expired
00247     if (prevTimeIndex_ == -1)
00248     {
00249         return false;
00250     }
00251     
00252     // force update
00253     prevTimeIndex_ = -1;
00254     return true;
00255 }
00256 
00257 
00258 bool Foam::sampledIsoSurfaceCell::update()
00259 {
00260     return updateGeometry();
00261 }
00262 
00263 
00264 Foam::tmp<Foam::scalarField>
00265 Foam::sampledIsoSurfaceCell::sample
00266 (
00267     const volScalarField& vField
00268 ) const
00269 {
00270     return sampleField(vField);
00271 }
00272 
00273 
00274 Foam::tmp<Foam::vectorField>
00275 Foam::sampledIsoSurfaceCell::sample
00276 (
00277     const volVectorField& vField
00278 ) const
00279 {
00280     return sampleField(vField);
00281 }
00282 
00283 
00284 Foam::tmp<Foam::sphericalTensorField>
00285 Foam::sampledIsoSurfaceCell::sample
00286 (
00287     const volSphericalTensorField& vField
00288 ) const
00289 {
00290     return sampleField(vField);
00291 }
00292 
00293 
00294 Foam::tmp<Foam::symmTensorField>
00295 Foam::sampledIsoSurfaceCell::sample
00296 (
00297     const volSymmTensorField& vField
00298 ) const
00299 {
00300     return sampleField(vField);
00301 }
00302 
00303 
00304 Foam::tmp<Foam::tensorField>
00305 Foam::sampledIsoSurfaceCell::sample
00306 (
00307     const volTensorField& vField
00308 ) const
00309 {
00310     return sampleField(vField);
00311 }
00312 
00313 
00314 Foam::tmp<Foam::scalarField>
00315 Foam::sampledIsoSurfaceCell::interpolate
00316 (
00317     const interpolation<scalar>& interpolator
00318 ) const
00319 {
00320     return interpolateField(interpolator);
00321 }
00322 
00323 
00324 Foam::tmp<Foam::vectorField>
00325 Foam::sampledIsoSurfaceCell::interpolate
00326 (
00327     const interpolation<vector>& interpolator
00328 ) const
00329 {
00330     return interpolateField(interpolator);
00331 }
00332 
00333 Foam::tmp<Foam::sphericalTensorField>
00334 Foam::sampledIsoSurfaceCell::interpolate
00335 (
00336     const interpolation<sphericalTensor>& interpolator
00337 ) const
00338 {
00339     return interpolateField(interpolator);
00340 }
00341 
00342 
00343 Foam::tmp<Foam::symmTensorField>
00344 Foam::sampledIsoSurfaceCell::interpolate
00345 (
00346     const interpolation<symmTensor>& interpolator
00347 ) const
00348 {
00349     return interpolateField(interpolator);
00350 }
00351 
00352 
00353 Foam::tmp<Foam::tensorField>
00354 Foam::sampledIsoSurfaceCell::interpolate
00355 (
00356     const interpolation<tensor>& interpolator
00357 ) const
00358 {
00359     return interpolateField(interpolator);
00360 }
00361 
00362 
00363 void Foam::sampledIsoSurfaceCell::print(Ostream& os) const
00364 {
00365     os  << "sampledIsoSurfaceCell: " << name() << " :"
00366         << "  field:" << isoField_
00367         << "  value:" << isoVal_;
00368         //<< "  faces:" << faces().size()   // possibly no geom yet
00369         //<< "  points:" << points().size();
00370 }
00371 
00372 
00373 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines