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

movingConeTopoFvMesh.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 "movingConeTopoFvMesh.H"
00027 #include <OpenFOAM/Time.H>
00028 #include <OpenFOAM/mapPolyMesh.H>
00029 #include <dynamicMesh/layerAdditionRemoval.H>
00030 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00031 #include <meshTools/meshTools.H>
00032 #include <OpenFOAM/OFstream.H>
00033 
00034 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00035 
00036 namespace Foam
00037 {
00038     defineTypeNameAndDebug(movingConeTopoFvMesh, 0);
00039 
00040     addToRunTimeSelectionTable
00041     (
00042         topoChangerFvMesh,
00043         movingConeTopoFvMesh,
00044         IOobject
00045     );
00046 }
00047 
00048 
00049 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00050 
00051 Foam::tmp<Foam::scalarField> Foam::movingConeTopoFvMesh::vertexMarkup
00052 (
00053     const pointField& p,
00054     const scalar& curLeft,
00055     const scalar& curRight
00056 ) const
00057 {
00058     Info<< "Updating vertex markup.  curLeft: "
00059         << curLeft << " curRight: " << curRight << endl;
00060 
00061     tmp<scalarField> tvertexMarkup(new scalarField(p.size()));
00062     scalarField& vertexMarkup = tvertexMarkup();
00063 
00064     forAll (p, pI)
00065     {
00066         if (p[pI].x() < curLeft - SMALL)
00067         {
00068             vertexMarkup[pI] = -1;
00069         }
00070         else if (p[pI].x() > curRight + SMALL)
00071         {
00072             vertexMarkup[pI] = 1;
00073         }
00074         else
00075         {
00076             vertexMarkup[pI] = 0;
00077         }
00078     }
00079 
00080     return tvertexMarkup;
00081 }
00082 
00083 
00084 void Foam::movingConeTopoFvMesh::addZonesAndModifiers()
00085 {
00086     // Add zones and modifiers for motion action
00087 
00088     if
00089     (
00090         pointZones().size()
00091      || faceZones().size()
00092      || cellZones().size()
00093      || topoChanger_.size()
00094     )
00095     {
00096         Info<< "void movingConeTopoFvMesh::addZonesAndModifiers() : "
00097             << "Zones and modifiers already present.  Skipping."
00098             << endl;
00099 
00100         return;
00101     }
00102 
00103     Info<< "Time = " << time().timeName() << endl
00104         << "Adding zones and modifiers to the mesh" << endl;
00105 
00106     const vectorField& fc = faceCentres();
00107     const vectorField& fa = faceAreas();
00108 
00109     labelList zone1(fc.size());
00110     boolList flipZone1(fc.size(), false);
00111     label nZoneFaces1 = 0;
00112 
00113     labelList zone2(fc.size());
00114     boolList flipZone2(fc.size(), false);
00115     label nZoneFaces2 = 0;
00116 
00117     forAll (fc, faceI)
00118     {
00119         if
00120         (
00121             fc[faceI].x() > -0.003501
00122          && fc[faceI].x() < -0.003499
00123         )
00124         {
00125             if ((fa[faceI] & vector(1, 0, 0)) < 0)
00126             {
00127                 flipZone1[nZoneFaces1] = true;
00128             }
00129 
00130             zone1[nZoneFaces1] = faceI;
00131             Info<< "face " << faceI << " for zone 1.  Flip: "
00132                 << flipZone1[nZoneFaces1] << endl;
00133             nZoneFaces1++;
00134         }
00135         else if
00136         (
00137             fc[faceI].x() > -0.00701
00138          && fc[faceI].x() < -0.00699
00139         )
00140         {
00141             zone2[nZoneFaces2] = faceI;
00142 
00143             if ((fa[faceI] & vector(1, 0, 0)) > 0)
00144             {
00145                 flipZone2[nZoneFaces2] = true;
00146             }
00147 
00148             Info << "face " << faceI << " for zone 2.  Flip: "
00149                 << flipZone2[nZoneFaces2] << endl;
00150             nZoneFaces2++;
00151         }
00152     }
00153 
00154     zone1.setSize(nZoneFaces1);
00155     flipZone1.setSize(nZoneFaces1);
00156 
00157     zone2.setSize(nZoneFaces2);
00158     flipZone2.setSize(nZoneFaces2);
00159 
00160     Info << "zone: " << zone1 << endl;
00161     Info << "zone: " << zone2 << endl;
00162 
00163     List<pointZone*> pz(0);
00164     List<faceZone*> fz(2);
00165     List<cellZone*> cz(0);
00166 
00167     label nFz = 0;
00168 
00169     fz[nFz] =
00170         new faceZone
00171         (
00172             "rightExtrusionFaces",
00173             zone1,
00174             flipZone1,
00175             nFz,
00176             faceZones()
00177         );
00178     nFz++;
00179 
00180     fz[nFz] =
00181         new faceZone
00182         (
00183             "leftExtrusionFaces",
00184             zone2,
00185             flipZone2,
00186             nFz,
00187             faceZones()
00188         );
00189     nFz++;
00190 
00191     fz.setSize(nFz);
00192 
00193     Info << "Adding mesh zones." << endl;
00194     addZones(pz, fz, cz);
00195 
00196 
00197     // Add layer addition/removal interfaces
00198 
00199     List<polyMeshModifier*> tm(2);
00200     label nMods = 0;
00201 
00202     tm[nMods] =
00203         new layerAdditionRemoval
00204         (
00205             "right",
00206             nMods,
00207             topoChanger_,
00208             "rightExtrusionFaces",
00209             readScalar
00210             (
00211                 motionDict_.subDict("right").lookup("minThickness")
00212             ),
00213             readScalar
00214             (
00215                 motionDict_.subDict("right").lookup("maxThickness")
00216             )
00217         );
00218     nMods++;
00219 
00220     tm[nMods] = new layerAdditionRemoval
00221     (
00222         "left",
00223         nMods,
00224         topoChanger_,
00225         "leftExtrusionFaces",
00226         readScalar
00227         (
00228             motionDict_.subDict("left").lookup("minThickness")
00229         ),
00230         readScalar
00231         (
00232             motionDict_.subDict("left").lookup("maxThickness")
00233         )
00234     );
00235     nMods++;
00236     tm.setSize(nMods);
00237 
00238     Info << "Adding " << nMods << " mesh modifiers" << endl;
00239     topoChanger_.addTopologyModifiers(tm);
00240 
00241     write();
00242 }
00243 
00244 
00245 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00246 
00247 // Construct from components
00248 Foam::movingConeTopoFvMesh::movingConeTopoFvMesh(const IOobject& io)
00249 :
00250     topoChangerFvMesh(io),
00251     motionDict_
00252     (
00253         IOdictionary
00254         (
00255             IOobject
00256             (
00257                 "dynamicMeshDict",
00258                 time().constant(),
00259                 *this,
00260                 IOobject::MUST_READ,
00261                 IOobject::NO_WRITE
00262             )
00263         ).subDict(typeName + "Coeffs")
00264     ),
00265     motionVelAmplitude_(motionDict_.lookup("motionVelAmplitude")),
00266     motionVelPeriod_(readScalar(motionDict_.lookup("motionVelPeriod"))),
00267     curMotionVel_
00268     (
00269         motionVelAmplitude_*
00270         Foam::sin(time().value()*M_PI/motionVelPeriod_)
00271     ),
00272     leftEdge_(readScalar(motionDict_.lookup("leftEdge"))),
00273     curLeft_(readScalar(motionDict_.lookup("leftObstacleEdge"))),
00274     curRight_(readScalar(motionDict_.lookup("rightObstacleEdge")))
00275 {
00276     Pout<< "Initial time:" << time().value()
00277         << " Initial curMotionVel_:" << curMotionVel_
00278         << endl;
00279 
00280     addZonesAndModifiers();
00281 
00282     curLeft_ = average
00283     (
00284         faceZones()
00285         [
00286             faceZones().findZoneID("leftExtrusionFaces")
00287         ]().localPoints()
00288     ).x() - SMALL;
00289 
00290     curRight_ = average
00291     (
00292         faceZones()
00293         [
00294             faceZones().findZoneID("rightExtrusionFaces")
00295         ]().localPoints()
00296     ).x() + SMALL;
00297 
00298     motionMask_ = vertexMarkup
00299     (
00300         points(),
00301         curLeft_,
00302         curRight_
00303     );
00304 }
00305 
00306 
00307 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00308 
00309 Foam::movingConeTopoFvMesh::~movingConeTopoFvMesh()
00310 {}
00311 
00312 
00313 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00314 
00315 bool Foam::movingConeTopoFvMesh::update()
00316 {
00317     // Do mesh changes (use inflation - put new points in topoChangeMap)
00318     autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
00319 
00320     // Calculate the new point positions depending on whether the
00321     // topological change has happened or not
00322     pointField newPoints;
00323 
00324     vector curMotionVel_ =
00325         motionVelAmplitude_*
00326         Foam::sin(time().value()*M_PI/motionVelPeriod_);
00327 
00328     Pout<< "time:" << time().value() << " curMotionVel_:" << curMotionVel_
00329         << " curLeft:" << curLeft_ << " curRight:" << curRight_
00330         << endl;
00331 
00332     if (topoChangeMap.valid())
00333     {
00334         Info << "Topology change. Calculating motion points" << endl;
00335 
00336         if (topoChangeMap().hasMotionPoints())
00337         {
00338             Info << "Topology change. Has premotion points" << endl;
00339             //Info<< "preMotionPoints:" << topoChangeMap().preMotionPoints()
00340             //    << endl;
00341 
00342             //mkDir(time().timePath());
00343             //{
00344             //    OFstream str(time().timePath()/"meshPoints.obj");
00345             //    Pout<< "Writing mesh with meshPoints to " << str.name()
00346             //        << endl;
00347             //
00348             //    const pointField& currentPoints = points();
00349             //    label vertI = 0;
00350             //    forAll(currentPoints, pointI)
00351             //    {
00352             //        meshTools::writeOBJ(str, currentPoints[pointI]);
00353             //        vertI++;
00354             //    }
00355             //    forAll(edges(), edgeI)
00356             //    {
00357             //        const edge& e = edges()[edgeI];
00358             //        str << "l " << e[0]+1 << ' ' << e[1]+1 << nl;
00359             //    }
00360             //}
00361             //{
00362             //    OFstream str(time().timePath()/"preMotionPoints.obj");
00363             //    Pout<< "Writing mesh with preMotionPoints to " << str.name()
00364             //        << endl;
00365             //
00366             //    const pointField& newPoints =
00367             //        topoChangeMap().preMotionPoints();
00368             //    label vertI = 0;
00369             //    forAll(newPoints, pointI)
00370             //    {
00371             //        meshTools::writeOBJ(str, newPoints[pointI]);
00372             //        vertI++;
00373             //    }
00374             //    forAll(edges(), edgeI)
00375             //    {
00376             //        const edge& e = edges()[edgeI];
00377             //        str << "l " << e[0]+1 << ' ' << e[1]+1 << nl;
00378             //    }
00379             //}
00380 
00381 
00382             motionMask_ =
00383                 vertexMarkup
00384                 (
00385                     topoChangeMap().preMotionPoints(),
00386                     curLeft_,
00387                     curRight_
00388                 );
00389 
00390             // Move points inside the motionMask
00391             newPoints =
00392                 topoChangeMap().preMotionPoints()
00393               + (
00394                     pos(0.5 - mag(motionMask_)) // cells above the body
00395                 )*curMotionVel_*time().deltaT().value();
00396         }
00397         else
00398         {
00399             Info << "Topology change. Already set mesh points" << endl;
00400 
00401             motionMask_ =
00402                 vertexMarkup
00403                 (
00404                     points(),
00405                     curLeft_,
00406                     curRight_
00407                 );
00408 
00409             // Move points inside the motionMask
00410             newPoints =
00411                 points()
00412               + (
00413                     pos(0.5 - mag(motionMask_)) // cells above the body
00414                 )*curMotionVel_*time().deltaT().value();
00415         }
00416     }
00417     else
00418     {
00419         Info << "No topology change" << endl;
00420         // Set the mesh motion
00421         newPoints =
00422             points()
00423           + (
00424                 pos(0.5 - mag(motionMask_)) // cells above the body
00425            )*curMotionVel_*time().deltaT().value();
00426     }
00427 
00428     // The mesh now contains the cells with zero volume
00429     Info << "Executing mesh motion" << endl;
00430     movePoints(newPoints);
00431     //  The mesh now has got non-zero volume cells
00432 
00433     curLeft_ = average
00434     (
00435         faceZones()
00436         [
00437             faceZones().findZoneID("leftExtrusionFaces")
00438         ]().localPoints()
00439     ).x() - SMALL;
00440 
00441     curRight_ = average
00442     (
00443         faceZones()
00444         [
00445             faceZones().findZoneID("rightExtrusionFaces")
00446         ]().localPoints()
00447     ).x() + SMALL;
00448 
00449 
00450     return true;
00451 }
00452 
00453 
00454 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines