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

distanceSurface.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 "distanceSurface.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 <sampling/isoSurface.H>
00033 // #include "isoSurfaceCell.H"
00034 
00035 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00036 
00037 namespace Foam
00038 {
00039     defineTypeNameAndDebug(distanceSurface, 0);
00040     addToRunTimeSelectionTable(sampledSurface, distanceSurface, word);
00041 }
00042 
00043 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00044 
00045 void Foam::distanceSurface::createGeometry()
00046 {
00047     if (debug)
00048     {
00049         Pout<< "distanceSurface::createGeometry :updating geometry." << endl;
00050     }
00051 
00052     // Clear any stored topologies
00053     facesPtr_.clear();
00054 
00055     // Clear derived data
00056     clearGeom();
00057 
00058     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
00059 
00060     // Distance to cell centres
00061     // ~~~~~~~~~~~~~~~~~~~~~~~~
00062 
00063     cellDistancePtr_.reset
00064     (
00065         new volScalarField
00066         (
00067             IOobject
00068             (
00069                 "cellDistance",
00070                 fvm.time().timeName(),
00071                 fvm.time(),
00072                 IOobject::NO_READ,
00073                 IOobject::NO_WRITE,
00074                 false
00075             ),
00076             fvm,
00077             dimensionedScalar("zero", dimLength, 0)
00078         )
00079     );
00080     volScalarField& cellDistance = cellDistancePtr_();
00081 
00082     // Internal field
00083     {
00084         const pointField& cc = fvm.C();
00085         scalarField& fld = cellDistance.internalField();
00086 
00087         List<pointIndexHit> nearest;
00088         surfPtr_().findNearest
00089         (
00090             cc,
00091             scalarField(cc.size(), GREAT),
00092             nearest
00093         );
00094 
00095         if (signed_)
00096         {
00097             vectorField normal;
00098             surfPtr_().getNormal(nearest, normal);
00099 
00100             forAll(nearest, i)
00101             {
00102                 vector d(cc[i]-nearest[i].hitPoint());
00103 
00104                 if ((d&normal[i]) > 0)
00105                 {
00106                     fld[i] = Foam::mag(d);
00107                 }
00108                 else
00109                 {
00110                     fld[i] = -Foam::mag(d);
00111                 }
00112             }
00113         }
00114         else
00115         {
00116             forAll(nearest, i)
00117             {
00118                 fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
00119             }
00120         }
00121     }
00122 
00123     // Patch fields
00124     {
00125         forAll(fvm.C().boundaryField(), patchI)
00126         {
00127             const pointField& cc = fvm.C().boundaryField()[patchI];
00128             fvPatchScalarField& fld = cellDistance.boundaryField()[patchI];
00129 
00130             List<pointIndexHit> nearest;
00131             surfPtr_().findNearest
00132             (
00133                 cc,
00134                 scalarField(cc.size(), GREAT),
00135                 nearest
00136             );
00137 
00138             if (signed_)
00139             {
00140                 vectorField normal;
00141                 surfPtr_().getNormal(nearest, normal);
00142 
00143                 forAll(nearest, i)
00144                 {
00145                     vector d(cc[i]-nearest[i].hitPoint());
00146 
00147                     if ((d&normal[i]) > 0)
00148                     {
00149                         fld[i] = Foam::mag(d);
00150                     }
00151                     else
00152                     {
00153                         fld[i] = -Foam::mag(d);
00154                     }
00155                 }
00156             }
00157             else
00158             {
00159                 forAll(nearest, i)
00160                 {
00161                     fld[i] = Foam::mag(cc[i] - nearest[i].hitPoint());
00162                 }
00163             }
00164         }
00165     }
00166 
00167 
00168     // On processor patches the mesh.C() will already be the cell centre
00169     // on the opposite side so no need to swap cellDistance.
00170 
00171 
00172     // Distance to points
00173     pointDistance_.setSize(fvm.nPoints());
00174     {
00175         const pointField& pts = fvm.points();
00176 
00177         List<pointIndexHit> nearest;
00178         surfPtr_().findNearest
00179         (
00180             pts,
00181             scalarField(pts.size(), GREAT),
00182             nearest
00183         );
00184 
00185         if (signed_)
00186         {
00187             vectorField normal;
00188             surfPtr_().getNormal(nearest, normal);
00189 
00190             forAll(nearest, i)
00191             {
00192                 vector d(pts[i]-nearest[i].hitPoint());
00193 
00194                 if ((d&normal[i]) > 0)
00195                 {
00196                     pointDistance_[i] = Foam::mag(d);
00197                 }
00198                 else
00199                 {
00200                     pointDistance_[i] = -Foam::mag(d);
00201                 }
00202             }
00203         }
00204         else
00205         {
00206             forAll(nearest, i)
00207             {
00208                 pointDistance_[i] = Foam::mag(pts[i]-nearest[i].hitPoint());
00209             }
00210         }
00211     }
00212 
00213 
00214     if (debug)
00215     {
00216         Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
00217         cellDistance.write();
00218         pointScalarField pDist
00219         (
00220             IOobject
00221             (
00222                 "pointDistance",
00223                 fvm.time().timeName(),
00224                 fvm.time(),
00225                 IOobject::NO_READ,
00226                 IOobject::NO_WRITE,
00227                 false
00228             ),
00229             pointMesh::New(fvm),
00230             dimensionedScalar("zero", dimLength, 0)
00231         );
00232         pDist.internalField() = pointDistance_;
00233 
00234         Pout<< "Writing point distance:" << pDist.objectPath() << endl;
00235         pDist.write();
00236     }
00237 
00238 
00239     //- Direct from cell field and point field.
00240     isoSurfPtr_.reset
00241     (
00242         new isoSurface
00243         (
00244             cellDistance,
00245             pointDistance_,
00246             distance_,
00247             regularise_
00248         )
00249     );
00250 
00251     if (debug)
00252     {
00253         print(Pout);
00254         Pout<< endl;
00255     }
00256 }
00257 
00258 
00259 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00260 
00261 Foam::distanceSurface::distanceSurface
00262 (
00263     const word& name,
00264     const polyMesh& mesh,
00265     const dictionary& dict
00266 )
00267 :
00268     sampledSurface(name, mesh, dict),
00269     surfPtr_
00270     (
00271         searchableSurface::New
00272         (
00273             dict.lookup("surfaceType"),
00274             IOobject
00275             (
00276                 dict.lookupOrDefault("surfaceName", name),  // name
00277                 mesh.time().constant(),                     // directory
00278                 "triSurface",                               // instance
00279                 mesh.time(),                                // registry
00280                 IOobject::MUST_READ,
00281                 IOobject::NO_WRITE
00282             ),
00283             dict
00284         )
00285     ),
00286     distance_(readScalar(dict.lookup("distance"))),
00287     signed_(readBool(dict.lookup("signed"))),
00288     regularise_(dict.lookupOrDefault("regularise", true)),
00289     zoneName_(word::null),
00290     needsUpdate_(true),
00291     isoSurfPtr_(NULL),
00292     facesPtr_(NULL)
00293 {
00294 //    dict.readIfPresent("zone", zoneName_);
00295 //
00296 //    if (debug && zoneName_.size())
00297 //    {
00298 //        if (mesh.cellZones().findZoneID(zoneName_) < 0)
00299 //        {
00300 //            Info<< "cellZone \"" << zoneName_
00301 //                << "\" not found - using entire mesh" << endl;
00302 //        }
00303 //    }
00304 }
00305 
00306 
00307 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00308 
00309 Foam::distanceSurface::~distanceSurface()
00310 {}
00311 
00312 
00313 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00314 
00315 bool Foam::distanceSurface::needsUpdate() const
00316 {
00317     return needsUpdate_;
00318 }
00319 
00320 
00321 bool Foam::distanceSurface::expire()
00322 {
00323     if (debug)
00324     {
00325         Pout<< "distanceSurface::expire :"
00326             << " have-facesPtr_:" << facesPtr_.valid()
00327             << " needsUpdate_:" << needsUpdate_ << endl;
00328     }
00329 
00330     // Clear any stored topologies
00331     facesPtr_.clear();
00332 
00333     // Clear derived data
00334     clearGeom();
00335 
00336     // already marked as expired
00337     if (needsUpdate_)
00338     {
00339         return false;
00340     }
00341 
00342     needsUpdate_ = true;
00343     return true;
00344 }
00345 
00346 
00347 bool Foam::distanceSurface::update()
00348 {
00349     if (debug)
00350     {
00351         Pout<< "distanceSurface::update :"
00352             << " have-facesPtr_:" << facesPtr_.valid()
00353             << " needsUpdate_:" << needsUpdate_ << endl;
00354     }
00355 
00356     if (!needsUpdate_)
00357     {
00358         return false;
00359     }
00360 
00361     createGeometry();
00362 
00363     needsUpdate_ = false;
00364     return true;
00365 }
00366 
00367 
00368 Foam::tmp<Foam::scalarField>
00369 Foam::distanceSurface::sample
00370 (
00371     const volScalarField& vField
00372 ) const
00373 {
00374     return sampleField(vField);
00375 }
00376 
00377 
00378 Foam::tmp<Foam::vectorField>
00379 Foam::distanceSurface::sample
00380 (
00381     const volVectorField& vField
00382 ) const
00383 {
00384     return sampleField(vField);
00385 }
00386 
00387 
00388 Foam::tmp<Foam::sphericalTensorField>
00389 Foam::distanceSurface::sample
00390 (
00391     const volSphericalTensorField& vField
00392 ) const
00393 {
00394     return sampleField(vField);
00395 }
00396 
00397 
00398 Foam::tmp<Foam::symmTensorField>
00399 Foam::distanceSurface::sample
00400 (
00401     const volSymmTensorField& vField
00402 ) const
00403 {
00404     return sampleField(vField);
00405 }
00406 
00407 
00408 Foam::tmp<Foam::tensorField>
00409 Foam::distanceSurface::sample
00410 (
00411     const volTensorField& vField
00412 ) const
00413 {
00414     return sampleField(vField);
00415 }
00416 
00417 
00418 Foam::tmp<Foam::scalarField>
00419 Foam::distanceSurface::interpolate
00420 (
00421     const interpolation<scalar>& interpolator
00422 ) const
00423 {
00424     return interpolateField(interpolator);
00425 }
00426 
00427 
00428 Foam::tmp<Foam::vectorField>
00429 Foam::distanceSurface::interpolate
00430 (
00431     const interpolation<vector>& interpolator
00432 ) const
00433 {
00434     return interpolateField(interpolator);
00435 }
00436 
00437 Foam::tmp<Foam::sphericalTensorField>
00438 Foam::distanceSurface::interpolate
00439 (
00440     const interpolation<sphericalTensor>& interpolator
00441 ) const
00442 {
00443     return interpolateField(interpolator);
00444 }
00445 
00446 
00447 Foam::tmp<Foam::symmTensorField>
00448 Foam::distanceSurface::interpolate
00449 (
00450     const interpolation<symmTensor>& interpolator
00451 ) const
00452 {
00453     return interpolateField(interpolator);
00454 }
00455 
00456 
00457 Foam::tmp<Foam::tensorField>
00458 Foam::distanceSurface::interpolate
00459 (
00460     const interpolation<tensor>& interpolator
00461 ) const
00462 {
00463     return interpolateField(interpolator);
00464 }
00465 
00466 
00467 void Foam::distanceSurface::print(Ostream& os) const
00468 {
00469     os  << "distanceSurface: " << name() << " :"
00470         << "  surface:" << surfPtr_().name()
00471         << "  distance:" << distance_
00472         << "  faces:" << faces().size()
00473         << "  points:" << points().size();
00474 }
00475 
00476 
00477 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines