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

mapLagrangian.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 "MapLagrangianFields.H"
00027 #include <lagrangian/Cloud.H>
00028 #include <lagrangian/passiveParticle.H>
00029 #include <meshTools/meshSearch.H>
00030 
00031 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00032 
00033 namespace Foam
00034 {
00035 
00036 static const scalar perturbFactor = 1E-6;
00037 
00038 
00039 // Special version of findCell that generates a cell guaranteed to be
00040 // compatible with tracking.
00041 static label findCell(const meshSearch& meshSearcher, const point& pt)
00042 {
00043     const polyMesh& mesh = meshSearcher.mesh();
00044 
00045     // Use tracking to find cell containing pt
00046     label cellI = meshSearcher.findCell(pt);
00047 
00048     if (cellI >= 0)
00049     {
00050         return cellI;
00051     }
00052     else
00053     {
00054         // See if particle on face by finding nearest face and shifting
00055         // particle.
00056 
00057         label faceI = meshSearcher.findNearestBoundaryFace(pt);
00058 
00059         if (faceI >= 0)
00060         {
00061             const point& cc = mesh.cellCentres()[mesh.faceOwner()[faceI]];
00062             const point perturbPt = (1-perturbFactor)*pt+perturbFactor*cc;
00063 
00064             return meshSearcher.findCell(perturbPt);
00065         }
00066     }
00067     return -1;
00068 }
00069 
00070 
00071 void mapLagrangian(const meshToMesh& meshToMeshInterp)
00072 {
00073     // Determine which particles are in meshTarget
00074     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00075 
00076     // target to source cell map
00077     const labelList& cellAddressing = meshToMeshInterp.cellAddressing();
00078 
00079     // Invert celladdressing to get source to target(s).
00080     // Note: could use sparse addressing but that is too storage inefficient
00081     // (Map<labelList>)
00082     labelListList sourceToTargets
00083     (
00084         invertOneToMany(meshToMeshInterp.fromMesh().nCells(), cellAddressing)
00085     );
00086 
00087     const fvMesh& meshSource = meshToMeshInterp.fromMesh();
00088     const fvMesh& meshTarget = meshToMeshInterp.toMesh();
00089     const pointField& targetCc = meshTarget.cellCentres();
00090 
00091 
00092     fileNameList cloudDirs
00093     (
00094         readDir
00095         (
00096             meshSource.time().timePath()/cloud::prefix,
00097             fileName::DIRECTORY
00098         )
00099     );
00100 
00101     forAll(cloudDirs, cloudI)
00102     {
00103         // Search for list of lagrangian objects for this time
00104         IOobjectList objects
00105         (
00106             meshSource,
00107             meshSource.time().timeName(),
00108             cloud::prefix/cloudDirs[cloudI]
00109         );
00110 
00111         IOobject* positionsPtr = objects.lookup("positions");
00112 
00113         if (positionsPtr)
00114         {
00115             Info<< nl << "    processing cloud " << cloudDirs[cloudI] << endl;
00116 
00117             // Read positions & cell
00118             Cloud<passiveParticle> sourceParcels
00119             (
00120                 meshSource,
00121                 cloudDirs[cloudI],
00122                 false
00123             );
00124             Info<< "    read " << sourceParcels.size()
00125                 << " parcels from source mesh." << endl;
00126 
00127             // Construct empty target cloud
00128             Cloud<passiveParticle> targetParcels
00129             (
00130                 meshTarget,
00131                 cloudDirs[cloudI],
00132                 IDLList<passiveParticle>()
00133             );
00134 
00135             label sourceParticleI = 0;
00136 
00137             // Indices of source particles that get added to targetParcels
00138             DynamicList<label> addParticles(sourceParcels.size());
00139 
00140             // Unmapped particles
00141             labelHashSet unmappedSource(sourceParcels.size());
00142 
00143 
00144             // Initial: track from fine-mesh cell centre to particle position
00145             // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00146             // This requires there to be no boundary in the way.
00147 
00148 
00149             forAllConstIter(Cloud<passiveParticle>, sourceParcels, iter)
00150             {
00151                 bool foundCell = false;
00152 
00153                 // Assume that cell from read parcel is the correct one...
00154                 if (iter().cell() >= 0)
00155                 {
00156                     const labelList& targetCells =
00157                         sourceToTargets[iter().cell()];
00158 
00159                     // Particle probably in one of the targetcells. Try
00160                     // all by tracking from their cell centre to the parcel
00161                     // position.
00162 
00163                     forAll(targetCells, i)
00164                     {
00165                         // Track from its cellcentre to position to make sure.
00166                         autoPtr<passiveParticle> newPtr
00167                         (
00168                             new passiveParticle
00169                             (
00170                                 targetParcels,
00171                                 targetCc[targetCells[i]],
00172                                 targetCells[i]
00173                             )
00174                         );
00175                         passiveParticle& newP = newPtr();
00176 
00177                         scalar fraction = 0;
00178                         label faceI = newP.track(iter().position(), fraction);
00179 
00180                         if (faceI < 0 && newP.cell() >= 0)
00181                         {
00182                             // Hit position.
00183                             foundCell = true;
00184                             addParticles.append(sourceParticleI);
00185                             targetParcels.addParticle(newPtr.ptr());
00186                             break;
00187                         }
00188                     }
00189                 }
00190 
00191                 if (!foundCell)
00192                 {
00193                     // Store for closer analysis
00194                     unmappedSource.insert(sourceParticleI);
00195                 }
00196 
00197                 sourceParticleI++;
00198             }
00199 
00200             Info<< "    after meshToMesh addressing found "
00201                 << targetParcels.size()
00202                 << " parcels in target mesh." << endl;
00203 
00204 
00205             // Do closer inspection for unmapped particles
00206             // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00207 
00208             if (unmappedSource.size())
00209             {
00210                 meshSearch targetSearcher(meshTarget, false);
00211 
00212                 sourceParticleI = 0;
00213 
00214                 forAllIter(Cloud<passiveParticle>, sourceParcels, iter)
00215                 {
00216                     if (unmappedSource.found(sourceParticleI))
00217                     {
00218                         label targetCell =
00219                             findCell(targetSearcher, iter().position());
00220 
00221                         if (targetCell >= 0)
00222                         {
00223                             unmappedSource.erase(sourceParticleI);
00224                             addParticles.append(sourceParticleI);
00225                 iter().cell()=targetCell;
00226                             targetParcels.addParticle
00227                             (
00228                                 sourceParcels.remove(&iter())
00229                             );
00230                         }
00231                     }
00232                     sourceParticleI++;
00233                 }
00234             }
00235             addParticles.shrink();
00236 
00237             Info<< "    after additional mesh searching found "
00238                 << targetParcels.size() << " parcels in target mesh." << endl;
00239 
00240             if (addParticles.size())
00241             {
00242                 IOPosition<passiveParticle>(targetParcels).write();
00243 
00244                 // addParticles now contains the indices of the sourceMesh
00245                 // particles that were appended to the target mesh.
00246 
00247                 // Map lagrangian fields
00248                 // ~~~~~~~~~~~~~~~~~~~~~
00249 
00250                 MapLagrangianFields<label>
00251                 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00252                 MapLagrangianFields<scalar>
00253                 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00254                 MapLagrangianFields<vector>
00255                 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00256                 MapLagrangianFields<sphericalTensor>
00257                 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00258                 MapLagrangianFields<symmTensor>
00259                 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00260                 MapLagrangianFields<tensor>
00261                 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00262             }
00263         }
00264     }
00265 }
00266 
00267 
00268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00269 
00270 } // End namespace Foam
00271 
00272 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines