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

fvMesh.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 "fvMesh.H"
00027 #include <finiteVolume/volFields.H>
00028 #include <finiteVolume/surfaceFields.H>
00029 #include <finiteVolume/slicedVolFields.H>
00030 #include <finiteVolume/slicedSurfaceFields.H>
00031 #include <OpenFOAM/SubField.H>
00032 #include <OpenFOAM/demandDrivenData.H>
00033 #include "fvMeshLduAddressing.H"
00034 #include <OpenFOAM/emptyPolyPatch.H>
00035 #include <OpenFOAM/mapPolyMesh.H>
00036 #include <finiteVolume/MapFvFields.H>
00037 #include <finiteVolume/fvMeshMapper.H>
00038 #include <OpenFOAM/mapClouds.H>
00039 
00040 #include <finiteVolume/volPointInterpolation.H>
00041 #include <finiteVolume/extendedLeastSquaresVectors.H>
00042 #include <finiteVolume/extendedLeastSquaresVectors.H>
00043 #include <finiteVolume/leastSquaresVectors.H>
00044 #include <finiteVolume/CentredFitData.H>
00045 #include <finiteVolume/linearFitPolynomial.H>
00046 #include <finiteVolume/quadraticFitPolynomial.H>
00047 #include <finiteVolume/quadraticLinearFitPolynomial.H>
00048 //#include "quadraticFitSnGradData.H"
00049 #include <finiteVolume/skewCorrectionVectors.H>
00050 
00051 
00052 #include <finiteVolume/centredCECCellToFaceStencilObject.H>
00053 #include <finiteVolume/centredCFCCellToFaceStencilObject.H>
00054 #include <finiteVolume/centredCPCCellToFaceStencilObject.H>
00055 #include <finiteVolume/centredFECCellToFaceStencilObject.H>
00056 #include <finiteVolume/upwindCECCellToFaceStencilObject.H>
00057 #include <finiteVolume/upwindCFCCellToFaceStencilObject.H>
00058 #include <finiteVolume/upwindCPCCellToFaceStencilObject.H>
00059 #include <finiteVolume/upwindFECCellToFaceStencilObject.H>
00060 
00061 #include <finiteVolume/centredCFCFaceToCellStencilObject.H>
00062 
00063 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00064 
00065 defineTypeNameAndDebug(Foam::fvMesh, 0);
00066 
00067 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00068 
00069 void Foam::fvMesh::clearGeomNotOldVol()
00070 {
00071     slicedVolScalarField::DimensionedInternalField* VPtr =
00072         static_cast<slicedVolScalarField::DimensionedInternalField*>(VPtr_);
00073     deleteDemandDrivenData(VPtr);
00074     VPtr_ = NULL;
00075 
00076     deleteDemandDrivenData(SfPtr_);
00077     deleteDemandDrivenData(magSfPtr_);
00078     deleteDemandDrivenData(CPtr_);
00079     deleteDemandDrivenData(CfPtr_);
00080 }
00081 
00082 
00083 void Foam::fvMesh::clearGeom()
00084 {
00085     clearGeomNotOldVol();
00086 
00087     deleteDemandDrivenData(V0Ptr_);
00088     deleteDemandDrivenData(V00Ptr_);
00089 
00090     // Mesh motion flux cannot be deleted here because the old-time flux
00091     // needs to be saved.
00092 
00093     // Things geometry dependent that are not updated.
00094     volPointInterpolation::Delete(*this);
00095     extendedLeastSquaresVectors::Delete(*this);
00096     leastSquaresVectors::Delete(*this);
00097     CentredFitData<linearFitPolynomial>::Delete(*this);
00098     CentredFitData<quadraticFitPolynomial>::Delete(*this);
00099     CentredFitData<quadraticLinearFitPolynomial>::Delete(*this);
00100     skewCorrectionVectors::Delete(*this);
00101     //quadraticFitSnGradData::Delete(*this);
00102 }
00103 
00104 
00105 void Foam::fvMesh::clearAddressing()
00106 {
00107     deleteDemandDrivenData(lduPtr_);
00108 
00109     // Hack until proper callbacks. Below are all the fvMesh-MeshObjects.
00110 
00111     volPointInterpolation::Delete(*this);
00112     extendedLeastSquaresVectors::Delete(*this);
00113     leastSquaresVectors::Delete(*this);
00114     CentredFitData<linearFitPolynomial>::Delete(*this);
00115     CentredFitData<quadraticFitPolynomial>::Delete(*this);
00116     CentredFitData<quadraticLinearFitPolynomial>::Delete(*this);
00117     skewCorrectionVectors::Delete(*this);
00118     //quadraticFitSnGradData::Delete(*this);
00119 
00120     centredCECCellToFaceStencilObject::Delete(*this);
00121     centredCFCCellToFaceStencilObject::Delete(*this);
00122     centredCPCCellToFaceStencilObject::Delete(*this);
00123     centredFECCellToFaceStencilObject::Delete(*this);
00124     // Is this geometry related - cells distorting to upwind direction?
00125     upwindCECCellToFaceStencilObject::Delete(*this);
00126     upwindCFCCellToFaceStencilObject::Delete(*this);
00127     upwindCPCCellToFaceStencilObject::Delete(*this);
00128     upwindFECCellToFaceStencilObject::Delete(*this);
00129 
00130     centredCFCFaceToCellStencilObject::Delete(*this);
00131 }
00132 
00133 
00134 void Foam::fvMesh::clearOut()
00135 {
00136     clearGeom();
00137     surfaceInterpolation::clearOut();
00138 
00139     clearAddressing();
00140 
00141     // Clear mesh motion flux
00142     deleteDemandDrivenData(phiPtr_);
00143 
00144     polyMesh::clearOut();
00145 }
00146 
00147 
00148 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00149 
00150 Foam::fvMesh::fvMesh(const IOobject& io)
00151 :
00152     polyMesh(io),
00153     surfaceInterpolation(*this),
00154     boundary_(*this, boundaryMesh()),
00155     lduPtr_(NULL),
00156     curTimeIndex_(time().timeIndex()),
00157     VPtr_(NULL),
00158     V0Ptr_(NULL),
00159     V00Ptr_(NULL),
00160     SfPtr_(NULL),
00161     magSfPtr_(NULL),
00162     CPtr_(NULL),
00163     CfPtr_(NULL),
00164     phiPtr_(NULL)
00165 {
00166     if (debug)
00167     {
00168         Info<< "Constructing fvMesh from IOobject"
00169             << endl;
00170     }
00171 
00172     // Check the existance of the cell volumes and read if present
00173     // and set the storage of V00
00174     if (isFile(time().timePath()/"V0"))
00175     {
00176         V0Ptr_ = new DimensionedField<scalar, volMesh>
00177         (
00178             IOobject
00179             (
00180                 "V0",
00181                 time().timeName(),
00182                 *this,
00183                 IOobject::MUST_READ,
00184                 IOobject::NO_WRITE
00185             ),
00186             *this
00187         );
00188 
00189         V00();
00190     }
00191 
00192     // Check the existance of the mesh fluxes, read if present and set the
00193     // mesh to be moving
00194     if (isFile(time().timePath()/"meshPhi"))
00195     {
00196         phiPtr_ = new surfaceScalarField
00197         (
00198             IOobject
00199             (
00200                 "meshPhi",
00201                 time().timeName(),
00202                 *this,
00203                 IOobject::MUST_READ,
00204                 IOobject::AUTO_WRITE
00205             ),
00206             *this
00207         );
00208 
00209         // The mesh is now considered moving so the old-time cell volumes
00210         // will be required for the time derivatives so if they haven't been
00211         // read initialise to the current cell volumes
00212         if (!V0Ptr_)
00213         {
00214             V0Ptr_ = new DimensionedField<scalar, volMesh>
00215             (
00216                 IOobject
00217                 (
00218                     "V0",
00219                     time().timeName(),
00220                     *this,
00221                     IOobject::NO_READ,
00222                     IOobject::NO_WRITE,
00223                     false
00224                 ),
00225                 V()
00226             );
00227         }
00228 
00229         moving(true);
00230     }
00231 }
00232 
00233 
00234 Foam::fvMesh::fvMesh
00235 (
00236     const IOobject& io,
00237     const Xfer<pointField>& points,
00238     const Xfer<faceList>& faces,
00239     const Xfer<labelList>& allOwner,
00240     const Xfer<labelList>& allNeighbour,
00241     const bool syncPar
00242 )
00243 :
00244     polyMesh(io, points, faces, allOwner, allNeighbour, syncPar),
00245     surfaceInterpolation(*this),
00246     boundary_(*this),
00247     lduPtr_(NULL),
00248     curTimeIndex_(time().timeIndex()),
00249     VPtr_(NULL),
00250     V0Ptr_(NULL),
00251     V00Ptr_(NULL),
00252     SfPtr_(NULL),
00253     magSfPtr_(NULL),
00254     CPtr_(NULL),
00255     CfPtr_(NULL),
00256     phiPtr_(NULL)
00257 {
00258     if (debug)
00259     {
00260         Info<< "Constructing fvMesh from components" << endl;
00261     }
00262 }
00263 
00264 
00265 Foam::fvMesh::fvMesh
00266 (
00267     const IOobject& io,
00268     const Xfer<pointField>& points,
00269     const Xfer<faceList>& faces,
00270     const Xfer<cellList>& cells,
00271     const bool syncPar
00272 )
00273 :
00274     polyMesh(io, points, faces, cells, syncPar),
00275     surfaceInterpolation(*this),
00276     boundary_(*this),
00277     lduPtr_(NULL),
00278     curTimeIndex_(time().timeIndex()),
00279     VPtr_(NULL),
00280     V0Ptr_(NULL),
00281     V00Ptr_(NULL),
00282     SfPtr_(NULL),
00283     magSfPtr_(NULL),
00284     CPtr_(NULL),
00285     CfPtr_(NULL),
00286     phiPtr_(NULL)
00287 {
00288     if (debug)
00289     {
00290         Info<< "Constructing fvMesh from components" << endl;
00291     }
00292 }
00293 
00294 
00295 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00296 
00297 Foam::fvMesh::~fvMesh()
00298 {
00299     clearOut();
00300 }
00301 
00302 
00303 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00304 
00305 void Foam::fvMesh::addFvPatches
00306 (
00307     const List<polyPatch*> & p,
00308     const bool validBoundary
00309 )
00310 {
00311     if (boundary().size())
00312     {
00313         FatalErrorIn
00314         (
00315             "fvMesh::addFvPatches(const List<polyPatch*>&, const bool)"
00316         )   << " boundary already exists"
00317             << abort(FatalError);
00318     }
00319 
00320     // first add polyPatches
00321     addPatches(p, validBoundary);
00322     boundary_.addPatches(boundaryMesh());
00323 }
00324 
00325 
00326 void Foam::fvMesh::removeFvBoundary()
00327 {
00328     if (debug)
00329     {
00330         Info<< "void fvMesh::removeFvBoundary(): "
00331             << "Removing boundary patches."
00332             << endl;
00333     }
00334 
00335     // Remove fvBoundaryMesh data first.
00336     boundary_.clear();
00337     boundary_.setSize(0);
00338     polyMesh::removeBoundary();
00339 
00340     clearOut();
00341 }
00342 
00343 
00344 Foam::polyMesh::readUpdateState Foam::fvMesh::readUpdate()
00345 {
00346     if (debug)
00347     {
00348         Info<< "polyMesh::readUpdateState fvMesh::readUpdate() : "
00349             << "Updating fvMesh.  ";
00350     }
00351 
00352     polyMesh::readUpdateState state = polyMesh::readUpdate();
00353 
00354     if (state == polyMesh::TOPO_PATCH_CHANGE)
00355     {
00356         if (debug)
00357         {
00358             Info << "Boundary and topological update" << endl;
00359         }
00360 
00361         boundary_.readUpdate(boundaryMesh());
00362 
00363         clearOut();
00364 
00365     }
00366     else if (state == polyMesh::TOPO_CHANGE)
00367     {
00368         if (debug)
00369         {
00370             Info << "Topological update" << endl;
00371         }
00372 
00373         clearOut();
00374     }
00375     else if (state == polyMesh::POINTS_MOVED)
00376     {
00377         if (debug)
00378         {
00379             Info << "Point motion update" << endl;
00380         }
00381 
00382         clearGeom();
00383     }
00384     else
00385     {
00386         if (debug)
00387         {
00388             Info << "No update" << endl;
00389         }
00390     }
00391 
00392     return state;
00393 }
00394 
00395 
00396 const Foam::fvBoundaryMesh& Foam::fvMesh::boundary() const
00397 {
00398     return boundary_;
00399 }
00400 
00401 
00402 const Foam::lduAddressing& Foam::fvMesh::lduAddr() const
00403 {
00404     if (!lduPtr_)
00405     {
00406         lduPtr_ = new fvMeshLduAddressing(*this);
00407     }
00408 
00409     return *lduPtr_;
00410 }
00411 
00412 
00413 void Foam::fvMesh::mapFields(const mapPolyMesh& meshMap)
00414 {
00415     // Create a mapper
00416     const fvMeshMapper mapper(*this, meshMap);
00417 
00418     // Map all the volFields in the objectRegistry
00419     MapGeometricFields<scalar, fvPatchField, fvMeshMapper, volMesh>
00420     (mapper);
00421     MapGeometricFields<vector, fvPatchField, fvMeshMapper, volMesh>
00422     (mapper);
00423     MapGeometricFields<sphericalTensor, fvPatchField, fvMeshMapper, volMesh>
00424     (mapper);
00425     MapGeometricFields<symmTensor, fvPatchField, fvMeshMapper, volMesh>
00426     (mapper);
00427     MapGeometricFields<tensor, fvPatchField, fvMeshMapper, volMesh>
00428     (mapper);
00429 
00430     // Map all the surfaceFields in the objectRegistry
00431     MapGeometricFields<scalar, fvsPatchField, fvMeshMapper, surfaceMesh>
00432     (mapper);
00433     MapGeometricFields<vector, fvsPatchField, fvMeshMapper, surfaceMesh>
00434     (mapper);
00435     MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
00436     (mapper);
00437     MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
00438     (mapper);
00439     MapGeometricFields<tensor, fvsPatchField, fvMeshMapper, surfaceMesh>
00440     (mapper);
00441 
00442     // Map all the clouds in the objectRegistry
00443     mapClouds(*this, meshMap);
00444 
00445 
00446     const labelList& cellMap = meshMap.cellMap();
00447 
00448     // Map the old volume. Just map to new cell labels.
00449     if (V0Ptr_)
00450     {
00451         scalarField& V0 = *V0Ptr_;
00452 
00453         scalarField savedV0(V0);
00454         V0.setSize(nCells());
00455 
00456         forAll(V0, i)
00457         {
00458             if (cellMap[i] > -1)
00459             {
00460                 V0[i] = savedV0[cellMap[i]];
00461             }
00462             else
00463             {
00464                 V0[i] = 0.0;
00465             }
00466         }
00467     }
00468 
00469     // Map the old-old volume. Just map to new cell labels.
00470     if (V00Ptr_)
00471     {
00472         scalarField& V00 = *V00Ptr_;
00473 
00474         scalarField savedV00(V00);
00475         V00.setSize(nCells());
00476 
00477         forAll(V00, i)
00478         {
00479             if (cellMap[i] > -1)
00480             {
00481                 V00[i] = savedV00[cellMap[i]];
00482             }
00483             else
00484             {
00485                 V00[i] = 0.0;
00486             }
00487         }
00488     }
00489 }
00490 
00491 
00492 // Temporary helper function to call move points on
00493 // MeshObjects
00494 template<class Type>
00495 void MeshObjectMovePoints(const Foam::fvMesh& mesh)
00496 {
00497     if (mesh.thisDb().foundObject<Type>(Type::typeName))
00498     {
00499         const_cast<Type&>
00500         (
00501             mesh.thisDb().lookupObject<Type>
00502             (
00503                 Type::typeName
00504             )
00505         ).movePoints();
00506     }
00507 }
00508 
00509 
00510 Foam::tmp<Foam::scalarField> Foam::fvMesh::movePoints(const pointField& p)
00511 {
00512     // Grab old time volumes if the time has been incremented
00513     if (curTimeIndex_ < time().timeIndex())
00514     {
00515         if (V00Ptr_ && V0Ptr_)
00516         {
00517             *V00Ptr_ = *V0Ptr_;
00518         }
00519 
00520         if (V0Ptr_)
00521         {
00522             *V0Ptr_ = V();
00523         }
00524         else
00525         {
00526             V0Ptr_ = new DimensionedField<scalar, volMesh>
00527             (
00528                 IOobject
00529                 (
00530                     "V0",
00531                     time().timeName(),
00532                     *this,
00533                     IOobject::NO_READ,
00534                     IOobject::NO_WRITE
00535                 ),
00536                 V()
00537             );
00538         }
00539 
00540         curTimeIndex_ = time().timeIndex();
00541     }
00542 
00543 
00544     // delete out of date geometrical information
00545     clearGeomNotOldVol();
00546 
00547 
00548     if (!phiPtr_)
00549     {
00550         // Create mesh motion flux
00551         phiPtr_ = new surfaceScalarField
00552         (
00553             IOobject
00554             (
00555                 "meshPhi",
00556                 this->time().timeName(),
00557                 *this,
00558                 IOobject::NO_READ,
00559                 IOobject::AUTO_WRITE
00560             ),
00561             *this,
00562             dimVolume/dimTime
00563         );
00564     }
00565     else
00566     {
00567         // Grab old time mesh motion fluxes if the time has been incremented
00568         if (phiPtr_->timeIndex() != time().timeIndex())
00569         {
00570             phiPtr_->oldTime();
00571         }
00572     }
00573 
00574     surfaceScalarField& phi = *phiPtr_;
00575 
00576     // Move the polyMesh and set the mesh motion fluxes to the swept-volumes
00577 
00578     scalar rDeltaT = 1.0/time().deltaT().value();
00579 
00580     tmp<scalarField> tsweptVols = polyMesh::movePoints(p);
00581     scalarField& sweptVols = tsweptVols();
00582 
00583     phi.internalField() = scalarField::subField(sweptVols, nInternalFaces());
00584     phi.internalField() *= rDeltaT;
00585 
00586     const fvPatchList& patches = boundary();
00587 
00588     forAll (patches, patchI)
00589     {
00590         phi.boundaryField()[patchI] = patches[patchI].patchSlice(sweptVols);
00591         phi.boundaryField()[patchI] *= rDeltaT;
00592     }
00593 
00594     boundary_.movePoints();
00595     surfaceInterpolation::movePoints();
00596 
00597 
00598     // Hack until proper callbacks. Below are all the fvMesh MeshObjects with a
00599     // movePoints function.
00600     MeshObjectMovePoints<volPointInterpolation>(*this);
00601     MeshObjectMovePoints<extendedLeastSquaresVectors>(*this);
00602     MeshObjectMovePoints<leastSquaresVectors>(*this);
00603     MeshObjectMovePoints<CentredFitData<linearFitPolynomial> >(*this);
00604     MeshObjectMovePoints<CentredFitData<quadraticFitPolynomial> >(*this);
00605     MeshObjectMovePoints<CentredFitData<quadraticLinearFitPolynomial> >(*this);
00606     MeshObjectMovePoints<skewCorrectionVectors>(*this);
00607     //MeshObjectMovePoints<quadraticFitSnGradData>(*this);
00608 
00609     return tsweptVols;
00610 }
00611 
00612 
00613 void Foam::fvMesh::updateMesh(const mapPolyMesh& mpm)
00614 {
00615     // Update polyMesh. This needs to keep volume existent!
00616     polyMesh::updateMesh(mpm);
00617 
00618     // Clear the sliced fields
00619     clearGeomNotOldVol();
00620 
00621     // Map all fields using current (i.e. not yet mapped) volume
00622     mapFields(mpm);
00623 
00624     // Clear the current volume and other geometry factors
00625     surfaceInterpolation::clearOut();
00626 
00627     clearAddressing();
00628 
00629     // handleMorph() should also clear out the surfaceInterpolation.
00630     // This is a temporary solution
00631     surfaceInterpolation::movePoints();
00632 }
00633 
00634 
00635 bool Foam::fvMesh::writeObjects
00636 (
00637     IOstream::streamFormat fmt,
00638     IOstream::versionNumber ver,
00639     IOstream::compressionType cmp
00640 ) const
00641 {
00642     return polyMesh::writeObject(fmt, ver, cmp);
00643 }
00644 
00645 
00646 //- Write mesh using IO settings from the time
00647 bool Foam::fvMesh::write() const
00648 {
00649     return polyMesh::write();
00650 }
00651 
00652 
00653 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
00654 
00655 bool Foam::fvMesh::operator!=(const fvMesh& bm) const
00656 {
00657     return &bm != this;
00658 }
00659 
00660 
00661 bool Foam::fvMesh::operator==(const fvMesh& bm) const
00662 {
00663     return &bm == this;
00664 }
00665 
00666 
00667 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines