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

sampledIsoSurface.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 "sampledIsoSurface.H"
00027 #include <OpenFOAM/dictionary.H>
00028 #include <finiteVolume/volFields.H>
00029 #include <finiteVolume/volPointInterpolation.H>
00030 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00031 
00032 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00033 
00034 namespace Foam
00035 {
00036     defineTypeNameAndDebug(sampledIsoSurface, 0);
00037     addNamedToRunTimeSelectionTable
00038     (
00039         sampledSurface,
00040         sampledIsoSurface,
00041         word,
00042         isoSurface
00043     );
00044 }
00045 
00046 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00047 
00048 void Foam::sampledIsoSurface::getIsoFields() const
00049 {
00050     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
00051 
00052     // Get volField
00053     // ~~~~~~~~~~~~
00054 
00055     if (fvm.foundObject<volScalarField>(isoField_))
00056     {
00057         if (debug)
00058         {
00059             Info<< "sampledIsoSurface::getIsoField() : lookup volField "
00060                 << isoField_ << endl;
00061         }
00062         storedVolFieldPtr_.clear();
00063         volFieldPtr_ = &fvm.lookupObject<volScalarField>(isoField_);
00064     }
00065     else
00066     {
00067         // Bit of a hack. Read field and store.
00068 
00069         if (debug)
00070         {
00071             Info<< "sampledIsoSurface::getIsoField() : checking "
00072                 << isoField_ << " for same time " << fvm.time().timeName()
00073                 << endl;
00074         }
00075 
00076         if
00077         (
00078             storedVolFieldPtr_.empty()
00079          || (fvm.time().timeName() != storedVolFieldPtr_().instance())
00080         )
00081         {
00082             if (debug)
00083             {
00084                 Info<< "sampledIsoSurface::getIsoField() : reading volField "
00085                     << isoField_ << " from time " << fvm.time().timeName()
00086                     << endl;
00087             }
00088 
00089             storedVolFieldPtr_.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             volFieldPtr_ = storedVolFieldPtr_.operator->();
00106         }
00107     }
00108 
00109 
00110 
00111     // Get pointField
00112     // ~~~~~~~~~~~~~~
00113 
00114     if (!subMeshPtr_.valid())
00115     {
00116         word pointFldName = "volPointInterpolate(" + isoField_ + ')';
00117 
00118         if (fvm.foundObject<pointScalarField>(pointFldName))
00119         {
00120             if (debug)
00121             {
00122                 Info<< "sampledIsoSurface::getIsoField() : lookup pointField "
00123                     << pointFldName << endl;
00124             }
00125             pointFieldPtr_ = &fvm.lookupObject<pointScalarField>(pointFldName);
00126         }
00127         else
00128         {
00129             // Not in registry. Interpolate.
00130 
00131             if (debug)
00132             {
00133                 Info<< "sampledIsoSurface::getIsoField() : checking pointField "
00134                     << pointFldName << " for same time "
00135                     << fvm.time().timeName() << endl;
00136             }
00137 
00138             if
00139             (
00140                 storedPointFieldPtr_.empty()
00141              || (fvm.time().timeName() != storedPointFieldPtr_().instance())
00142             )
00143             {
00144                 if (debug)
00145                 {
00146                     Info<< "sampledIsoSurface::getIsoField() :"
00147                         << " interpolating volField " << volFieldPtr_->name()
00148                         << " to get pointField " << pointFldName << endl;
00149                 }
00150 
00151                 storedPointFieldPtr_.reset
00152                 (
00153                     volPointInterpolation::New(fvm)
00154                     .interpolate(*volFieldPtr_).ptr()
00155                 );
00156                 storedPointFieldPtr_->checkOut();
00157                 pointFieldPtr_ = storedPointFieldPtr_.operator->();
00158             }
00159         }
00160 
00161 
00162         // If averaging redo the volField. Can only be done now since needs the
00163         // point field.
00164         if (average_)
00165         {
00166             storedVolFieldPtr_.reset(average(fvm, *pointFieldPtr_).ptr());
00167             volFieldPtr_ = storedVolFieldPtr_.operator->();
00168         }
00169 
00170 
00171         if (debug)
00172         {
00173             Info<< "sampledIsoSurface::getIsoField() : volField "
00174                 << volFieldPtr_->name() << " min:" << min(*volFieldPtr_).value()
00175                 << " max:" << max(*volFieldPtr_).value() << endl;
00176             Info<< "sampledIsoSurface::getIsoField() : pointField "
00177                 << pointFieldPtr_->name()
00178                 << " min:" << gMin(pointFieldPtr_->internalField())
00179                 << " max:" << gMax(pointFieldPtr_->internalField()) << endl;
00180         }
00181     }
00182     else
00183     {
00184         // Get subMesh variants
00185         const fvMesh& subFvm = subMeshPtr_().subMesh();
00186 
00187         // Either lookup on the submesh or subset the whole-mesh volField
00188 
00189         if (subFvm.foundObject<volScalarField>(isoField_))
00190         {
00191             if (debug)
00192             {
00193                 Info<< "sampledIsoSurface::getIsoField() :"
00194                     << " submesh lookup volField "
00195                     << isoField_ << endl;
00196             }
00197             storedVolSubFieldPtr_.clear();
00198             volSubFieldPtr_ = &subFvm.lookupObject<volScalarField>(isoField_);
00199         }
00200         else
00201         {
00202             if (debug)
00203             {
00204                 Info<< "sampledIsoSurface::getIsoField() : subsetting volField "
00205                     << isoField_ << endl;
00206             }
00207             storedVolSubFieldPtr_.reset
00208             (
00209                 subMeshPtr_().interpolate
00210                 (
00211                     *volFieldPtr_
00212                 ).ptr()
00213             );
00214             storedVolSubFieldPtr_->checkOut();
00215             volSubFieldPtr_ = storedVolSubFieldPtr_.operator->();
00216         }
00217 
00218 
00219         // Pointfield on submesh
00220 
00221         word pointFldName =
00222             "volPointInterpolate("
00223           + volSubFieldPtr_->name()
00224           + ')';
00225 
00226         if (subFvm.foundObject<pointScalarField>(pointFldName))
00227         {
00228             if (debug)
00229             {
00230                 Info<< "sampledIsoSurface::getIsoField() :"
00231                     << " submesh lookup pointField " << pointFldName << endl;
00232             }
00233             storedPointSubFieldPtr_.clear();
00234             pointSubFieldPtr_ = &subFvm.lookupObject<pointScalarField>
00235             (
00236                 pointFldName
00237             );
00238         }
00239         else
00240         {
00241             if (debug)
00242             {
00243                 Info<< "sampledIsoSurface::getIsoField() :"
00244                     << " interpolating submesh volField "
00245                     << volSubFieldPtr_->name()
00246                     << " to get submesh pointField " << pointFldName << endl;
00247             }
00248             storedPointSubFieldPtr_.reset
00249             (
00250                 volPointInterpolation::New
00251                 (
00252                     subFvm
00253                 ).interpolate(*volSubFieldPtr_).ptr()
00254             );
00255             storedPointSubFieldPtr_->checkOut();
00256             pointSubFieldPtr_ = storedPointSubFieldPtr_.operator->();
00257         }
00258 
00259 
00260 
00261         // If averaging redo the volField. Can only be done now since needs the
00262         // point field.
00263         if (average_)
00264         {
00265             storedVolSubFieldPtr_.reset
00266             (
00267                 average(subFvm, *pointSubFieldPtr_).ptr()
00268             );
00269             volSubFieldPtr_ = storedVolSubFieldPtr_.operator->();
00270         }
00271 
00272 
00273         if (debug)
00274         {
00275             Info<< "sampledIsoSurface::getIsoField() : volSubField "
00276                 << volSubFieldPtr_->name()
00277                 << " min:" << min(*volSubFieldPtr_).value()
00278                 << " max:" << max(*volSubFieldPtr_).value() << endl;
00279             Info<< "sampledIsoSurface::getIsoField() : pointSubField "
00280                 << pointSubFieldPtr_->name()
00281                 << " min:" << gMin(pointSubFieldPtr_->internalField())
00282                 << " max:" << gMax(pointSubFieldPtr_->internalField()) << endl;
00283         }
00284     }
00285 }
00286 
00287 
00288 Foam::tmp<Foam::volScalarField> Foam::sampledIsoSurface::average
00289 (
00290     const fvMesh& mesh,
00291     const pointScalarField& pfld
00292 ) const
00293 {
00294     tmp<volScalarField> tcellAvg
00295     (
00296         new volScalarField
00297         (
00298             IOobject
00299             (
00300                 "cellAvg",
00301                 mesh.time().timeName(),
00302                 mesh,
00303                 IOobject::NO_READ,
00304                 IOobject::NO_WRITE,
00305                 false
00306             ),
00307             mesh,
00308             dimensionedScalar("zero", dimless, scalar(0.0))
00309         )
00310     );
00311     volScalarField& cellAvg = tcellAvg();
00312 
00313     labelField nPointCells(mesh.nCells(), 0);
00314     {
00315         for (label pointI = 0; pointI < mesh.nPoints(); pointI++)
00316         {
00317             const labelList& pCells = mesh.pointCells(pointI);
00318 
00319             forAll(pCells, i)
00320             {
00321                 label cellI = pCells[i];
00322 
00323                 cellAvg[cellI] += pfld[pointI];
00324                 nPointCells[cellI]++;
00325             }
00326         }
00327     }
00328     forAll(cellAvg, cellI)
00329     {
00330         cellAvg[cellI] /= nPointCells[cellI];
00331     }
00332     // Give value to calculatedFvPatchFields
00333     cellAvg.correctBoundaryConditions();
00334 
00335     return tcellAvg;
00336 }
00337 
00338 
00339 Foam::tmp<Foam::pointScalarField> Foam::sampledIsoSurface::average
00340 (
00341     const pointMesh& pMesh,
00342     const volScalarField& fld
00343 ) const
00344 {
00345     tmp<pointScalarField> tpointAvg
00346     (
00347         new pointScalarField
00348         (
00349             IOobject
00350             (
00351                 "pointAvg",
00352                 fld.time().timeName(),
00353                 fld.db(),
00354                 IOobject::NO_READ,
00355                 IOobject::NO_WRITE,
00356                 false
00357             ),
00358             pMesh,
00359             dimensionedScalar("zero", dimless, scalar(0.0))
00360         )
00361     );
00362     pointScalarField& pointAvg = tpointAvg();
00363 
00364     for (label pointI = 0; pointI < fld.mesh().nPoints(); pointI++)
00365     {
00366         const labelList& pCells = fld.mesh().pointCells(pointI);
00367 
00368         forAll(pCells, i)
00369         {
00370             pointAvg[pointI] += fld[pCells[i]];
00371         }
00372         pointAvg[pointI] /= pCells.size();
00373     }
00374     // Give value to calculatedFvPatchFields
00375     pointAvg.correctBoundaryConditions();
00376 
00377     return tpointAvg;
00378 }
00379 
00380 
00381 bool Foam::sampledIsoSurface::updateGeometry() const
00382 {
00383     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
00384 
00385     // no update needed
00386     if (fvm.time().timeIndex() == prevTimeIndex_)
00387     {
00388         return false;
00389     }
00390 
00391     // Get any subMesh
00392     if (zoneID_.index() != -1 && !subMeshPtr_.valid())
00393     {
00394         const polyBoundaryMesh& patches = mesh().boundaryMesh();
00395 
00396         // Patch to put exposed internal faces into
00397         label exposedPatchI = patches.findPatchID(exposedPatchName_);
00398 
00399         if (debug)
00400         {
00401             Info<< "Allocating subset of size "
00402                 << mesh().cellZones()[zoneID_.index()].size()
00403                 << " with exposed faces into patch "
00404                 << patches[exposedPatchI].name() << endl;
00405         }
00406 
00407         subMeshPtr_.reset
00408         (
00409             new fvMeshSubset(fvm)
00410         );
00411         subMeshPtr_().setLargeCellSubset
00412         (
00413             labelHashSet(mesh().cellZones()[zoneID_.index()]),
00414             exposedPatchI
00415         );
00416     }
00417 
00418 
00419     prevTimeIndex_ = fvm.time().timeIndex();
00420     getIsoFields();
00421 
00422     // Clear any stored topo
00423     surfPtr_.clear();
00424     facesPtr_.clear();
00425 
00426     // Clear derived data
00427     clearGeom();
00428 
00429     if (subMeshPtr_.valid())
00430     {
00431         surfPtr_.reset
00432         (
00433             new isoSurface
00434             (
00435                 *volSubFieldPtr_,
00436                 *pointSubFieldPtr_,
00437                 isoVal_,
00438                 regularise_,
00439                 mergeTol_
00440             )
00441         );
00442     }
00443     else
00444     {
00445         surfPtr_.reset
00446         (
00447             new isoSurface
00448             (
00449                 *volFieldPtr_,
00450                 *pointFieldPtr_,
00451                 isoVal_,
00452                 regularise_,
00453                 mergeTol_
00454             )
00455         );
00456     }
00457 
00458 
00459     if (debug)
00460     {
00461         Pout<< "sampledIsoSurface::updateGeometry() : constructed iso:"
00462             << nl
00463             << "    regularise     : " << regularise_ << nl
00464             << "    average        : " << average_ << nl
00465             << "    isoField       : " << isoField_ << nl
00466             << "    isoValue       : " << isoVal_ << nl;
00467         if (subMeshPtr_.valid())
00468         {
00469             Pout<< "    zone size      : " << subMeshPtr_().subMesh().nCells()
00470                 << nl;
00471         }
00472         Pout<< "    points         : " << points().size() << nl
00473             << "    tris           : " << surface().size() << nl
00474             << "    cut cells      : " << surface().meshCells().size()
00475             << endl;
00476     }
00477 
00478     return true;
00479 }
00480 
00481 
00482 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00483 
00484 Foam::sampledIsoSurface::sampledIsoSurface
00485 (
00486     const word& name,
00487     const polyMesh& mesh,
00488     const dictionary& dict
00489 )
00490 :
00491     sampledSurface(name, mesh, dict),
00492     isoField_(dict.lookup("isoField")),
00493     isoVal_(readScalar(dict.lookup("isoValue"))),
00494     mergeTol_(dict.lookupOrDefault("mergeTol", 1E-6)),
00495     regularise_(dict.lookupOrDefault("regularise", true)),
00496     average_(dict.lookupOrDefault("average", false)),
00497     zoneID_(dict.lookupOrDefault("zone", word::null), mesh.cellZones()),
00498     exposedPatchName_(word::null),
00499     surfPtr_(NULL),
00500     facesPtr_(NULL),
00501     prevTimeIndex_(-1),
00502     storedVolFieldPtr_(NULL),
00503     volFieldPtr_(NULL),
00504     storedPointFieldPtr_(NULL),
00505     pointFieldPtr_(NULL)
00506 {
00507     if (!sampledSurface::interpolate())
00508     {
00509         FatalIOErrorIn
00510         (
00511             "sampledIsoSurface::sampledIsoSurface"
00512             "(const word&, const polyMesh&, const dictionary&)",
00513             dict
00514         )   << "Non-interpolated iso surface not supported since triangles"
00515             << " span across cells." << exit(FatalIOError);
00516     }
00517 
00518     if (zoneID_.index() != -1)
00519     {
00520         dict.lookup("exposedPatchName") >> exposedPatchName_;
00521 
00522         if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
00523         {
00524             FatalIOErrorIn
00525             (
00526                 "sampledIsoSurface::sampledIsoSurface"
00527                 "(const word&, const polyMesh&, const dictionary&)",
00528                 dict
00529             )   << "Cannot find patch " << exposedPatchName_
00530                 << " in which to put exposed faces." << endl
00531                 << "Valid patches are " << mesh.boundaryMesh().names()
00532                 << exit(FatalIOError);
00533         }
00534 
00535         if (debug && zoneID_.index() != -1)
00536         {
00537             Info<< "Restricting to cellZone " << zoneID_.name()
00538                 << " with exposed internal faces into patch "
00539                 << exposedPatchName_ << endl;
00540         }
00541     }
00542 }
00543 
00544 
00545 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00546 
00547 Foam::sampledIsoSurface::~sampledIsoSurface()
00548 {}
00549 
00550 
00551 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00552 
00553 bool Foam::sampledIsoSurface::needsUpdate() const
00554 {
00555     const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
00556 
00557     return fvm.time().timeIndex() != prevTimeIndex_;
00558 }
00559 
00560 
00561 bool Foam::sampledIsoSurface::expire()
00562 {
00563     surfPtr_.clear();
00564     facesPtr_.clear();
00565     subMeshPtr_.clear();
00566 
00567     // Clear derived data
00568     clearGeom();
00569 
00570     // already marked as expired
00571     if (prevTimeIndex_ == -1)
00572     {
00573         return false;
00574     }
00575 
00576     // force update
00577     prevTimeIndex_ = -1;
00578     return true;
00579 }
00580 
00581 
00582 bool Foam::sampledIsoSurface::update()
00583 {
00584     return updateGeometry();
00585 }
00586 
00587 
00588 Foam::tmp<Foam::scalarField> Foam::sampledIsoSurface::sample
00589 (
00590     const volScalarField& vField
00591 ) const
00592 {
00593     return sampleField(vField);
00594 }
00595 
00596 
00597 Foam::tmp<Foam::vectorField> Foam::sampledIsoSurface::sample
00598 (
00599     const volVectorField& vField
00600 ) const
00601 {
00602     return sampleField(vField);
00603 }
00604 
00605 
00606 Foam::tmp<Foam::sphericalTensorField> Foam::sampledIsoSurface::sample
00607 (
00608     const volSphericalTensorField& vField
00609 ) const
00610 {
00611     return sampleField(vField);
00612 }
00613 
00614 
00615 Foam::tmp<Foam::symmTensorField> Foam::sampledIsoSurface::sample
00616 (
00617     const volSymmTensorField& vField
00618 ) const
00619 {
00620     return sampleField(vField);
00621 }
00622 
00623 
00624 Foam::tmp<Foam::tensorField> Foam::sampledIsoSurface::sample
00625 (
00626     const volTensorField& vField
00627 ) const
00628 {
00629     return sampleField(vField);
00630 }
00631 
00632 
00633 Foam::tmp<Foam::scalarField> Foam::sampledIsoSurface::interpolate
00634 (
00635     const interpolation<scalar>& interpolator
00636 ) const
00637 {
00638     return interpolateField(interpolator);
00639 }
00640 
00641 
00642 Foam::tmp<Foam::vectorField> Foam::sampledIsoSurface::interpolate
00643 (
00644     const interpolation<vector>& interpolator
00645 ) const
00646 {
00647     return interpolateField(interpolator);
00648 }
00649 
00650 Foam::tmp<Foam::sphericalTensorField> Foam::sampledIsoSurface::interpolate
00651 (
00652     const interpolation<sphericalTensor>& interpolator
00653 ) const
00654 {
00655     return interpolateField(interpolator);
00656 }
00657 
00658 
00659 Foam::tmp<Foam::symmTensorField> Foam::sampledIsoSurface::interpolate
00660 (
00661     const interpolation<symmTensor>& interpolator
00662 ) const
00663 {
00664     return interpolateField(interpolator);
00665 }
00666 
00667 
00668 Foam::tmp<Foam::tensorField> Foam::sampledIsoSurface::interpolate
00669 (
00670     const interpolation<tensor>& interpolator
00671 ) const
00672 {
00673     return interpolateField(interpolator);
00674 }
00675 
00676 
00677 void Foam::sampledIsoSurface::print(Ostream& os) const
00678 {
00679     os  << "sampledIsoSurface: " << name() << " :"
00680         << "  field   :" << isoField_
00681         << "  value   :" << isoVal_;
00682         //<< "  faces:" << faces().size()       // note: possibly no geom yet
00683         //<< "  points:" << points().size();
00684 }
00685 
00686 
00687 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines