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

interpolationCellPointFace.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 Description
00025 
00026 \*---------------------------------------------------------------------------*/
00027 
00028 #include "interpolationCellPointFace.H"
00029 #include <finiteVolume/volFields.H>
00030 #include <OpenFOAM/polyMesh.H>
00031 #include <finiteVolume/volPointInterpolation.H>
00032 #include <finiteVolume/linear.H>
00033 #include "findCellPointFaceTet.H"
00034 #include "findCellPointFaceTriangle.H"
00035 
00036 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00037 
00038 namespace Foam
00039 {
00040 
00041 // * * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * //
00042 
00043 template<class Type>
00044 interpolationCellPointFace<Type>::interpolationCellPointFace
00045 (
00046     const GeometricField<Type, fvPatchField, volMesh>& psi
00047 )
00048 :
00049     interpolation<Type>(psi),
00050     psip_(volPointInterpolation::New(psi.mesh()).interpolate(psi)),
00051     psis_(linearInterpolate(psi))
00052 {}
00053 
00054 
00055 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00056 
00057 template<class Type>
00058 Type interpolationCellPointFace<Type>::interpolate
00059 (
00060     const vector& position,
00061     const label nCell,
00062     const label facei
00063 ) const
00064 {
00065     Type ts[4];
00066     vector tetPoints[4];
00067     scalar phi[4], phiCandidate[4];
00068     label tetLabelCandidate[2], tetPointLabels[2];
00069 
00070     Type t = pTraits<Type>::zero;
00071 
00072     // only use face information when the position is on a face
00073     if (facei < 0)
00074     {
00075         const vector& cellCentre = this->pMesh_.cellCentres()[nCell];
00076         const labelList& cellFaces = this->pMesh_.cells()[nCell];
00077         
00078         vector projection = position - cellCentre;
00079         tetPoints[3] = cellCentre;
00080         
00081         // *********************************************************************
00082         // project the cell-center through the point onto the face
00083         // and get the closest face, ...
00084         // *********************************************************************
00085         
00086         bool foundTet = false;
00087         label closestFace = -1;
00088         scalar minDistance = GREAT;
00089         
00090         forAll(cellFaces, facei)
00091         {
00092             label nFace = cellFaces[facei];
00093         
00094             vector normal = this->pMeshFaceAreas_[nFace];
00095             normal /= mag(normal);
00096         
00097             const vector& faceCentreTmp = this->pMeshFaceCentres_[nFace];
00098         
00099             scalar multiplierNumerator = (faceCentreTmp - cellCentre) & normal;
00100             scalar multiplierDenominator = projection & normal;
00101         
00102             // if normal and projection are not orthogonal this could be the one...
00103             if (mag(multiplierDenominator) > SMALL)
00104             {
00105                 scalar multiplier = multiplierNumerator/multiplierDenominator;
00106                 vector iPoint = cellCentre + multiplier*projection;
00107                 scalar dist = mag(position - iPoint);
00108         
00109                 if (dist < minDistance)
00110                 {
00111                     closestFace = nFace;
00112                     minDistance = dist;
00113                 }
00114             }
00115         }
00116         
00117         // *********************************************************************
00118         // find the tetrahedron containing 'position'
00119         // from the cell center, face center and
00120         // two other points on the face
00121         // *********************************************************************
00122         
00123         minDistance = GREAT;
00124         if (closestFace != -1)
00125         {
00126             label nFace = closestFace;
00127             foundTet = findTet
00128             (
00129                 position,
00130                 nFace,
00131                 tetPoints,
00132                 tetLabelCandidate,
00133                 tetPointLabels,
00134                 phi,
00135                 phiCandidate,
00136                 closestFace,
00137                 minDistance
00138             );
00139         }
00140         
00141         if (!foundTet)
00142         {
00143             // check if the position is 'just' outside a tet
00144             if (minDistance < 1.0e-5)
00145             {
00146                 foundTet = true;
00147                 for (label i=0; i<4; i++)
00148                 {
00149                     phi[i] = phiCandidate[i];
00150                 }
00151                 tetPointLabels[0] = tetLabelCandidate[0];
00152                 tetPointLabels[1] = tetLabelCandidate[1];
00153             }
00154         }
00155         
00156         // *********************************************************************
00157         // if the search failed check all the cell-faces
00158         // *********************************************************************
00159         
00160         if (!foundTet)
00161         {
00162             minDistance = GREAT;
00163         
00164             label facei = 0;
00165             while (facei < cellFaces.size() && !foundTet)
00166             {
00167                 label nFace = cellFaces[facei];
00168                 if (nFace < this->pMeshFaceAreas_.size())
00169                 {
00170                     foundTet = findTet
00171                     (
00172                         position,
00173                         nFace,
00174                         tetPoints,
00175                         tetLabelCandidate,
00176                         tetPointLabels,
00177                         phi,
00178                         phiCandidate,
00179                         closestFace,
00180                         minDistance
00181                     );
00182                 }
00183                 facei++;
00184             }
00185         }
00186         
00187         if (!foundTet)
00188         {
00189             // check if the position is 'just' outside a tet
00190             // this time with a more tolerant limit
00191             if (minDistance < 1.0e-3)
00192             {
00193                 foundTet = true;
00194                 for (label i=0; i<4; i++)
00195                 {
00196                     phi[i] = phiCandidate[i];
00197                 }
00198                 tetPointLabels[0] = tetLabelCandidate[0];
00199                 tetPointLabels[1] = tetLabelCandidate[1];
00200             }
00201         }
00202         
00203         // *********************************************************************
00204         // if the tet was found do the interpolation,
00205         // otherwise use the closest face value
00206         // *********************************************************************
00207         
00208         if (foundTet)
00209         {
00210             for (label i=0; i<2; i++)
00211             {
00212                 ts[i] = psip_[tetPointLabels[i]];
00213             }
00214         
00215             if (closestFace < psis_.size())
00216             {
00217                 ts[2] = psis_[closestFace];
00218             }
00219             else
00220             {
00221                 label patchi = this->pMesh_.boundaryMesh().whichPatch(closestFace);
00222         
00223                 // If the boundary patch is not empty use the face value
00224                 // else use the cell value
00225                 if (this->psi_.boundaryField()[patchi].size())
00226                 {
00227                     ts[2] = this->psi_.boundaryField()[patchi]
00228                         [this->pMesh_.boundaryMesh()[patchi].whichFace(closestFace)];
00229                 }
00230                 else
00231                 {
00232                     ts[2] = this->psi_[nCell];
00233                 }
00234             }
00235             
00236             ts[3] = this->psi_[nCell];
00237             
00238             for (label n=0; n<4; n++)
00239             {
00240                 phi[n] = min(1.0, phi[n]);
00241                 phi[n] = max(0.0, phi[n]);
00242         
00243                 t += phi[n]*ts[n];
00244             }
00245         }
00246         else
00247         {
00248             Info<< "interpolationCellPointFace<Type>::interpolate"
00249                 << "(const vector&, const label nCell) const : "
00250                 << "search failed; using closest cellFace value" << endl
00251                 << "cell number " << nCell << tab
00252                 << "position " << position << endl;
00253         
00254             if (closestFace < psis_.size())
00255             {
00256                 t = psis_[closestFace];
00257             }
00258             else
00259             {
00260                 label patchi = this->pMesh_.boundaryMesh().whichPatch(closestFace);
00261         
00262                 // If the boundary patch is not empty use the face value
00263                 // else use the cell value
00264                 if (this->psi_.boundaryField()[patchi].size())
00265                 {
00266                     t = this->psi_.boundaryField()[patchi]
00267                         [this->pMesh_.boundaryMesh()[patchi].whichFace(closestFace)];
00268                 }
00269                 else
00270                 {
00271                     t = this->psi_[nCell];
00272                 }
00273             }
00274         }
00275     }
00276     else
00277     {
00278         bool foundTriangle = findTriangle
00279         (
00280             position,
00281             facei,
00282             tetPointLabels,
00283             phi
00284         );
00285 
00286         if (foundTriangle)
00287         {
00288             // add up the point values ...
00289             for (label i=0; i<2; i++)
00290             {
00291                 Type vel = psip_[tetPointLabels[i]];
00292                 t += phi[i]*vel;
00293             }
00294 
00295             // ... and the face value
00296             if (facei < psis_.size())
00297             {
00298                 t += phi[2]*psis_[facei];
00299             }
00300             else
00301             {
00302                 label patchi = this->pMesh_.boundaryMesh().whichPatch(facei);
00303                 
00304                 // If the boundary patch is not empty use the face value
00305                 // else use the cell value
00306                 if (this->psi_.boundaryField()[patchi].size())
00307                 {
00308                     t += phi[2]*this->psi_.boundaryField()[patchi]
00309                         [this->pMesh_.boundaryMesh()[patchi].whichFace(facei)];
00310                 }
00311                 else
00312                 {
00313                     t += phi[2]*this->psi_[nCell];
00314                 }
00315             }
00316         }
00317         else
00318         {
00319             // use face value only
00320             if (facei < psis_.size())
00321             {
00322                 t = psis_[facei];
00323             }
00324             else
00325             {
00326                 label patchi = this->pMesh_.boundaryMesh().whichPatch(facei);
00327                 
00328                 // If the boundary patch is not empty use the face value
00329                 // else use the cell value
00330                 if (this->psi_.boundaryField()[patchi].size())
00331                 {
00332                     t = this->psi_.boundaryField()[patchi]
00333                         [this->pMesh_.boundaryMesh()[patchi].whichFace(facei)];
00334                 }
00335                 else
00336                 {
00337                     t = this->psi_[nCell];
00338                 }
00339             }
00340         }
00341     }
00342 
00343     return t;
00344 }
00345 
00346 
00347 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00348 
00349 } // End namespace Foam
00350 
00351 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines