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

calculateMeshToMeshAddressing.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     private member of meshToMesh.
00026     Calculates mesh to mesh addressing pattern (for each cell from one mesh
00027     find the closest cell centre in the other mesh).
00028 
00029 \*---------------------------------------------------------------------------*/
00030 
00031 #include "meshToMesh.H"
00032 #include <OpenFOAM/SubField.H>
00033 
00034 #include <meshTools/octree.H>
00035 #include <meshTools/octreeDataCell.H>
00036 #include <meshTools/octreeDataFace.H>
00037 
00038 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00039 
00040 namespace Foam
00041 {
00042 
00043 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00044 
00045 void meshToMesh::calcAddressing()
00046 {
00047     if (debug)
00048     {
00049         Info<< "meshToMesh::calculateAddressing() : "
00050             << "calculating mesh-to-mesh cell addressing" << endl;
00051     }
00052 
00053     // set reference to cells
00054     const cellList& fromCells = fromMesh_.cells();
00055     const pointField& fromPoints = fromMesh_.points();
00056 
00057     // In an attempt to preserve the efficiency of linear search and prevent
00058     // failure, a RESCUE mechanism will be set up. Here, we shall mark all
00059     // cells next to the solid boundaries. If such a cell is found as the
00060     // closest, the relationship between the origin and cell will be examined.
00061     // If the origin is outside the cell, a global n-squared search is
00062     // triggered.
00063 
00064     // SETTING UP RESCUE
00065 
00066     // visit all boundaries and mark the cell next to the boundary.
00067 
00068     if (debug)
00069     {
00070         Info<< "meshToMesh::calculateAddressing() : "
00071             << "Setting up rescue" << endl;
00072     }
00073 
00074     List<bool> boundaryCell(fromCells.size(), false);
00075 
00076     // set reference to boundary
00077     const polyPatchList& patchesFrom = fromMesh_.boundaryMesh();
00078 
00079     forAll (patchesFrom, patchI)
00080     {
00081         // get reference to cells next to the boundary
00082         const unallocLabelList& bCells = patchesFrom[patchI].faceCells();
00083 
00084         forAll (bCells, faceI)
00085         {
00086             boundaryCell[bCells[faceI]] = true;
00087         }
00088     }
00089 
00090     treeBoundBox meshBb(fromPoints);
00091 
00092     scalar typDim = meshBb.avgDim()/(2.0*cbrt(scalar(fromCells.size())));
00093 
00094     treeBoundBox shiftedBb
00095     (
00096         meshBb.min(),
00097         meshBb.max() + vector(typDim, typDim, typDim)
00098     );
00099 
00100     if (debug)
00101     {
00102         Info<< "\nMesh" << endl;
00103         Info<< "   bounding box           : " << meshBb << endl;
00104         Info<< "   bounding box (shifted) : " << shiftedBb << endl;
00105         Info<< "   typical dimension      :" << shiftedBb.typDim() << endl;
00106     }
00107 
00108     // Wrap indices and mesh information into helper object
00109     octreeDataCell shapes(fromMesh_);
00110 
00111     octree<octreeDataCell> oc
00112     (
00113         shiftedBb,  // overall bounding box
00114         shapes,     // all information needed to do checks on cells
00115         1,          // min. levels
00116         10.0,       // max. size of leaves
00117         10.0        // maximum ratio of cubes v.s. cells
00118     );
00119 
00120     if (debug)
00121     {
00122         oc.printStats(Info);
00123     }
00124 
00125     cellAddresses
00126     (
00127         cellAddressing_,
00128         toMesh_.cellCentres(),
00129         fromMesh_,
00130         boundaryCell,
00131         oc
00132     );
00133 
00134     forAll (toMesh_.boundaryMesh(), patchi)
00135     {
00136         const polyPatch& toPatch = toMesh_.boundaryMesh()[patchi];
00137 
00138         if (cuttingPatches_.found(toPatch.name()))
00139         {
00140             boundaryAddressing_[patchi].setSize(toPatch.size());
00141 
00142             cellAddresses
00143             (
00144                 boundaryAddressing_[patchi],
00145                 toPatch.faceCentres(),
00146                 fromMesh_,
00147                 boundaryCell,
00148                 oc
00149             );
00150         }
00151         else if
00152         (
00153             patchMap_.found(toPatch.name())
00154          && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
00155         )
00156         {
00157             const polyPatch& fromPatch = fromMesh_.boundaryMesh()
00158             [
00159                 fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
00160             ];
00161 
00162             if (fromPatch.empty())
00163             {
00164                 WarningIn("meshToMesh::calcAddressing()")
00165                     << "Source patch " << fromPatch.name()
00166                     << " has no faces. Not performing mapping for it."
00167                     << endl;
00168                 boundaryAddressing_[patchi] = -1;
00169             }
00170             else
00171             {
00172                 treeBoundBox wallBb(fromPatch.localPoints());
00173                 scalar typDim =
00174                     wallBb.avgDim()/(2.0*sqrt(scalar(fromPatch.size())));
00175 
00176                 treeBoundBox shiftedBb
00177                 (
00178                     wallBb.min(),
00179                     wallBb.max() + vector(typDim, typDim, typDim)
00180                 );
00181 
00182                 // Wrap data for octree into container
00183                 octreeDataFace shapes(fromPatch);
00184 
00185                 octree<octreeDataFace> oc
00186                 (
00187                     shiftedBb,  // overall search domain
00188                     shapes,     // all information needed to do checks on cells
00189                     1,          // min levels
00190                     20.0,       // maximum ratio of cubes v.s. cells
00191                     2.0
00192                 );
00193 
00194 
00195                 const vectorField::subField centresToBoundary =
00196                     toPatch.faceCentres();
00197 
00198                 boundaryAddressing_[patchi].setSize(toPatch.size());
00199 
00200                 treeBoundBox tightest;
00201                 scalar tightestDist;
00202 
00203                 forAll(toPatch, toi)
00204                 {
00205                     tightest = wallBb;                  // starting search bb
00206                     tightestDist = Foam::GREAT;        // starting max distance
00207 
00208                     boundaryAddressing_[patchi][toi] = oc.findNearest
00209                     (
00210                         centresToBoundary[toi],
00211                         tightest,
00212                         tightestDist
00213                     );
00214                 }
00215             }
00216         }
00217     }
00218 
00219     if (debug)
00220     {
00221         Info<< "meshToMesh::calculateAddressing() : "
00222             << "finished calculating mesh-to-mesh acell ddressing" << endl;
00223     }
00224 }
00225 
00226 
00227 void meshToMesh::cellAddresses
00228 (
00229     labelList& cellAddressing_,
00230     const pointField& points,
00231     const fvMesh& fromMesh,
00232     const List<bool>& boundaryCell,
00233     const octree<octreeDataCell>& oc
00234 ) const
00235 {
00236     // the implemented search method is a simple neighbour array search.
00237     // It starts from a cell zero, searches its neighbours and finds one
00238     // which is nearer to the target point than the current position.
00239     // The location of the "current position" is reset to that cell and
00240     // search through the neighbours continues. The search is finished
00241     // when all the neighbours of the cell are farther from the target
00242     // point than the current cell
00243 
00244     // set curCell label to zero (start)
00245     register label curCell = 0;
00246 
00247     // set reference to cell to cell addressing
00248     const vectorField& centresFrom = fromMesh.cellCentres();
00249     const labelListList& cc = fromMesh.cellCells();
00250 
00251     forAll (points, toI)
00252     {
00253         // pick up target position
00254         const vector& p = points[toI];
00255 
00256         // set the sqr-distance
00257         scalar distSqr = magSqr(p - centresFrom[curCell]);
00258 
00259         bool closer;
00260 
00261         do
00262         {
00263             closer = false;
00264 
00265             // set the current list of neighbouring cells
00266             const labelList& neighbours = cc[curCell];
00267 
00268             forAll (neighbours, nI)
00269             {
00270                 scalar curDistSqr =
00271                     magSqr(p - centresFrom[neighbours[nI]]);
00272 
00273                 // search through all the neighbours.
00274                 // If the cell is closer, reset current cell and distance
00275                 if (curDistSqr < (1 - SMALL)*distSqr)
00276                 {
00277                     curCell = neighbours[nI];
00278                     distSqr = curDistSqr;
00279                     closer = true;    // a closer neighbour has been found
00280                 }
00281             }
00282         } while (closer);
00283 
00284         cellAddressing_[toI] = -1;
00285 
00286         // Check point is actually in the nearest cell
00287         if (fromMesh.pointInCell(p, curCell))
00288         {
00289             cellAddressing_[toI] = curCell;
00290         }
00291         else
00292         {
00293             // If curCell is a boundary cell then the point maybe either outside
00294             // the domain or in an other region of the doamin, either way use
00295             // the octree search to find it.
00296             if (boundaryCell[curCell])
00297             {
00298                 cellAddressing_[toI] = oc.find(p);
00299             }
00300             else
00301             {
00302                 // If not on the boundary search the neighbours
00303                 bool found = false;
00304 
00305                 // set the current list of neighbouring cells
00306                 const labelList& neighbours = cc[curCell];
00307 
00308                 forAll (neighbours, nI)
00309                 {
00310                     // search through all the neighbours.
00311                     // If point is in neighbour reset current cell
00312                     if (fromMesh.pointInCell(p, neighbours[nI]))
00313                     {
00314                         cellAddressing_[toI] = neighbours[nI];
00315                         found = true;
00316                         break;
00317                     }
00318                 }
00319 
00320                 if (!found)
00321                 {
00322                     // If still not found search the neighbour-neighbours
00323 
00324                     // set the current list of neighbouring cells
00325                     const labelList& neighbours = cc[curCell];
00326 
00327                     forAll (neighbours, nI)
00328                     {
00329                         // set the current list of neighbour-neighbouring cells
00330                         const labelList& nn = cc[neighbours[nI]];
00331 
00332                         forAll (nn, nI)
00333                         {
00334                             // search through all the neighbours.
00335                             // If point is in neighbour reset current cell
00336                             if (fromMesh.pointInCell(p, nn[nI]))
00337                             {
00338                                 cellAddressing_[toI] = nn[nI];
00339                                 found = true;
00340                                 break;
00341                             }
00342                         }
00343                         if (found) break;
00344                     }
00345                 }
00346 
00347                 if (!found)
00348                 {
00349                     // Still not found so us the octree
00350                     cellAddressing_[toI] = oc.find(p);
00351                 }
00352             }
00353         }
00354     }
00355 }
00356 
00357 
00358 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00359 
00360 } // End namespace Foam
00361 
00362 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines