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

mergeOrSplitBaffles.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 Application
00025     mergeOrSplitBaffles
00026 
00027 Description
00028     Detects faces that share points (baffles). Either merge them or
00029     duplicate the points.
00030 
00031 Note
00032     - can only handle pairwise boundary faces. So three faces using
00033       the same points is not handled (is illegal mesh anyway)
00034 
00035     - there is no option to only split/merge some baffles.
00036 
00037     - surfaces consisting of duplicate faces can be topologically split
00038     if the points on the interior of the surface cannot walk to all the
00039     cells that use them in one go.
00040 
00041     - Parallel operation (where duplicate face is perpendicular to a coupled
00042     boundary) is supported but not really tested.
00043     (Note that coupled faces themselves are not seen as duplicate faces)
00044 
00045 Usage
00046 
00047     - mergeOrSplitBaffles [OPTIONS]
00048 
00049     @param -split \n
00050     Split duplicate surfaces.
00051 
00052     @param -overwrite \n
00053     Overwrite existing data.
00054 
00055     @param -detectOnly \n
00056     Do no processing.
00057 
00058     @param -case <dir>\n
00059     Case directory.
00060 
00061     @param -parallel \n
00062     Run in parallel.
00063 
00064     @param -help \n
00065     Display help message.
00066 
00067     @param -doc \n
00068     Display Doxygen API documentation page for this application.
00069 
00070     @param -srcDoc \n
00071     Display Doxygen source documentation page for this application.
00072 
00073 \*---------------------------------------------------------------------------*/
00074 
00075 #include <OpenFOAM/argList.H>
00076 #include <OpenFOAM/Time.H>
00077 #include <OpenFOAM/syncTools.H>
00078 #include <meshTools/faceSet.H>
00079 #include <meshTools/pointSet.H>
00080 #include <meshTools/meshTools.H>
00081 #include <dynamicMesh/polyTopoChange.H>
00082 #include <dynamicMesh/polyRemoveFace.H>
00083 #include <dynamicMesh/polyModifyFace.H>
00084 #include <OpenFOAM/indirectPrimitivePatch.H>
00085 #include <OpenFOAM/processorPolyPatch.H>
00086 #include <dynamicMesh/localPointRegion.H>
00087 #include <dynamicMesh/duplicatePoints.H>
00088 #include <OpenFOAM/ReadFields.H>
00089 #include <finiteVolume/volFields.H>
00090 #include <finiteVolume/surfaceFields.H>
00091 
00092 using namespace Foam;
00093 
00094 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00095 
00096 void insertDuplicateMerge
00097 (
00098     const polyMesh& mesh,
00099     const labelList& duplicates,
00100     polyTopoChange& meshMod
00101 )
00102 {
00103     const faceList& faces = mesh.faces();
00104     const labelList& faceOwner = mesh.faceOwner();
00105     const faceZoneMesh& faceZones = mesh.faceZones();
00106 
00107     forAll(duplicates, bFaceI)
00108     {
00109         label otherFaceI = duplicates[bFaceI];
00110 
00111         if (otherFaceI != -1 && otherFaceI > bFaceI)
00112         {
00113             // Two duplicate faces. Merge.
00114 
00115             label face0 = mesh.nInternalFaces() + bFaceI;
00116             label face1 = mesh.nInternalFaces() + otherFaceI;
00117 
00118             label own0 = faceOwner[face0];
00119             label own1 = faceOwner[face1];
00120 
00121             if (own0 < own1)
00122             {
00123                 // Use face0 as the new internal face.
00124                 label zoneID = faceZones.whichZone(face0);
00125                 bool zoneFlip = false;
00126 
00127                 if (zoneID >= 0)
00128                 {
00129                     const faceZone& fZone = faceZones[zoneID];
00130                     zoneFlip = fZone.flipMap()[fZone.whichFace(face0)];
00131                 }
00132 
00133                 meshMod.setAction(polyRemoveFace(face1));
00134                 meshMod.setAction
00135                 (
00136                     polyModifyFace
00137                     (
00138                         faces[face0],           // modified face
00139                         face0,                  // label of face being modified
00140                         own0,                   // owner
00141                         own1,                   // neighbour
00142                         false,                  // face flip
00143                         -1,                     // patch for face
00144                         false,                  // remove from zone
00145                         zoneID,                 // zone for face
00146                         zoneFlip                // face flip in zone
00147                     )
00148                 );
00149             }
00150             else
00151             {
00152                 // Use face1 as the new internal face.
00153                 label zoneID = faceZones.whichZone(face1);
00154                 bool zoneFlip = false;
00155 
00156                 if (zoneID >= 0)
00157                 {
00158                     const faceZone& fZone = faceZones[zoneID];
00159                     zoneFlip = fZone.flipMap()[fZone.whichFace(face1)];
00160                 }
00161 
00162                 meshMod.setAction(polyRemoveFace(face0));
00163                 meshMod.setAction
00164                 (
00165                     polyModifyFace
00166                     (
00167                         faces[face1],           // modified face
00168                         face1,                  // label of face being modified
00169                         own1,                   // owner
00170                         own0,                   // neighbour
00171                         false,                  // face flip
00172                         -1,                     // patch for face
00173                         false,                  // remove from zone
00174                         zoneID,                 // zone for face
00175                         zoneFlip                // face flip in zone
00176                     )
00177                 );
00178             }
00179         }
00180     }
00181 }
00182 
00183 
00184 labelList findBaffles(const polyMesh& mesh, const labelList& boundaryFaces)
00185 {
00186     // Get all duplicate face labels (in boundaryFaces indices!).
00187     labelList duplicates = localPointRegion::findDuplicateFaces
00188     (
00189         mesh,
00190         boundaryFaces
00191     );
00192 
00193 
00194     // Check that none are on processor patches
00195     const polyBoundaryMesh& patches = mesh.boundaryMesh();
00196 
00197     forAll(duplicates, bFaceI)
00198     {
00199         if (duplicates[bFaceI] != -1)
00200         {
00201             label faceI = mesh.nInternalFaces() + bFaceI;
00202             label patchI = patches.whichPatch(faceI);
00203 
00204             if (isA<processorPolyPatch>(patches[patchI]))
00205             {
00206                 FatalErrorIn("findBaffles(const polyMesh&, const labelList&)")
00207                     << "Duplicate face " << faceI
00208                     << " is on a processorPolyPatch."
00209                     << "This is not allowed." << nl
00210                     << "Face:" << faceI
00211                     << " is on patch:" << patches[patchI].name()
00212                     << abort(FatalError);
00213             }
00214         }
00215     }
00216 
00217 
00218     // Write to faceSet for ease of postprocessing.
00219     {
00220         faceSet duplicateSet
00221         (
00222             mesh,
00223             "duplicateFaces",
00224             (mesh.nFaces() - mesh.nInternalFaces())/256
00225         );
00226 
00227         forAll(duplicates, bFaceI)
00228         {
00229             label otherFaceI = duplicates[bFaceI];
00230 
00231             if (otherFaceI != -1 && otherFaceI > bFaceI)
00232             {
00233                 duplicateSet.insert(mesh.nInternalFaces() + bFaceI);
00234                 duplicateSet.insert(mesh.nInternalFaces() + otherFaceI);
00235             }
00236         }
00237 
00238         Pout<< "Writing " << duplicateSet.size()
00239             << " duplicate faces to faceSet " << duplicateSet.objectPath()
00240             << nl << endl;
00241         duplicateSet.write();
00242     }
00243 
00244     return duplicates;
00245 }
00246 
00247 
00248 
00249 
00250 int main(int argc, char *argv[])
00251 {
00252 #   include <OpenFOAM/addRegionOption.H>
00253     argList::validOptions.insert("split", "");
00254     argList::validOptions.insert("overwrite", "");
00255     argList::validOptions.insert("detectOnly", "");
00256 #   include <OpenFOAM/setRootCase.H>
00257 #   include <OpenFOAM/createTime.H>
00258     runTime.functionObjects().off();
00259 #   include <OpenFOAM/createNamedMesh.H>
00260     const word oldInstance = mesh.pointsInstance();
00261 
00262     bool split      = args.optionFound("split");
00263     bool overwrite  = args.optionFound("overwrite");
00264     bool detectOnly = args.optionFound("detectOnly");
00265 
00266     // Collect all boundary faces
00267     labelList boundaryFaces(mesh.nFaces() - mesh.nInternalFaces());
00268 
00269     forAll(boundaryFaces, i)
00270     {
00271         boundaryFaces[i] = i+mesh.nInternalFaces();
00272     }
00273 
00274 
00275     if (detectOnly)
00276     {
00277         findBaffles(mesh, boundaryFaces);
00278 
00279         return 0;
00280     }
00281 
00282 
00283 
00284     // Read objects in time directory
00285     IOobjectList objects(mesh, runTime.timeName());
00286 
00287     // Read vol fields.
00288 
00289     PtrList<volScalarField> vsFlds;
00290     ReadFields(mesh, objects, vsFlds);
00291 
00292     PtrList<volVectorField> vvFlds;
00293     ReadFields(mesh, objects, vvFlds);
00294 
00295     PtrList<volSphericalTensorField> vstFlds;
00296     ReadFields(mesh, objects, vstFlds);
00297 
00298     PtrList<volSymmTensorField> vsymtFlds;
00299     ReadFields(mesh, objects, vsymtFlds);
00300 
00301     PtrList<volTensorField> vtFlds;
00302     ReadFields(mesh, objects, vtFlds);
00303 
00304     // Read surface fields.
00305 
00306     PtrList<surfaceScalarField> ssFlds;
00307     ReadFields(mesh, objects, ssFlds);
00308 
00309     PtrList<surfaceVectorField> svFlds;
00310     ReadFields(mesh, objects, svFlds);
00311 
00312     PtrList<surfaceSphericalTensorField> sstFlds;
00313     ReadFields(mesh, objects, sstFlds);
00314 
00315     PtrList<surfaceSymmTensorField> ssymtFlds;
00316     ReadFields(mesh, objects, ssymtFlds);
00317 
00318     PtrList<surfaceTensorField> stFlds;
00319     ReadFields(mesh, objects, stFlds);
00320 
00321 
00322     // Mesh change engine
00323     polyTopoChange meshMod(mesh);
00324 
00325 
00326     if (split)
00327     {
00328         Pout<< "Topologically splitting duplicate surfaces"
00329             << ", i.e. duplicating points internal to duplicate surfaces."
00330             << nl << endl;
00331 
00332         // Analyse which points need to be duplicated
00333         localPointRegion regionSide(mesh);
00334 
00335         // Point duplication engine
00336         duplicatePoints pointDuplicator(mesh);
00337 
00338         // Insert topo changes
00339         pointDuplicator.setRefinement(regionSide, meshMod);
00340     }
00341     else
00342     {
00343         Pout<< "Merging duplicate faces."
00344             << nl << endl;
00345 
00346         // Get all duplicate face labels (in boundaryFaces indices!).
00347         labelList duplicates(findBaffles(mesh, boundaryFaces));
00348 
00349         // Merge into internal faces.
00350         insertDuplicateMerge(mesh, duplicates, meshMod);
00351     }
00352 
00353     if (!overwrite)
00354     {
00355         runTime++;
00356     }
00357 
00358     // Change the mesh. No inflation.
00359     autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
00360 
00361     // Update fields
00362     mesh.updateMesh(map);
00363 
00364     // Move mesh (since morphing does not do this)
00365     if (map().hasMotionPoints())
00366     {
00367         mesh.movePoints(map().preMotionPoints());
00368     }
00369 
00370     if (overwrite)
00371     {
00372         mesh.setInstance(oldInstance);
00373     }
00374     Pout<< "Writing mesh to time " << runTime.timeName() << endl;
00375     mesh.write();
00376 
00377     // Dump duplicated points (if any)
00378     if (split)
00379     {
00380         const labelList& pointMap = map().pointMap();
00381 
00382         labelList nDupPerPoint(map().nOldPoints(), 0);
00383 
00384         pointSet dupPoints(mesh, "duplicatedPoints", 100);
00385 
00386         forAll(pointMap, pointI)
00387         {
00388             label oldPointI = pointMap[pointI];
00389 
00390             nDupPerPoint[oldPointI]++;
00391 
00392             if (nDupPerPoint[oldPointI] > 1)
00393             {
00394                 dupPoints.insert(map().reversePointMap()[oldPointI]);
00395                 dupPoints.insert(pointI);
00396             }
00397         }
00398 
00399         Pout<< "Writing " << dupPoints.size()
00400             << " duplicated points to pointSet "
00401             << dupPoints.objectPath() << nl << endl;
00402 
00403         dupPoints.write();
00404     }
00405 
00406     Info<< "End\n" << endl;
00407 
00408     return 0;
00409 }
00410 
00411 
00412 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines