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

smoothDelta.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 "smoothDelta.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <OpenFOAM/FaceCellWave.H>
00029 
00030 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00031 
00032 namespace Foam
00033 {
00034 
00035 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00036 
00037 defineTypeNameAndDebug(smoothDelta, 0);
00038 addToRunTimeSelectionTable(LESdelta, smoothDelta, dictionary);
00039 
00040 scalar smoothDelta::deltaData::maxDeltaRatio = 1.2;
00041 
00042 
00043 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00044 
00045 // Fill changedFaces (with face labels) and changedFacesInfo (with delta)
00046 // This is the initial set of faces from which to start the waves.
00047 // Since there might be lots of places with delta jumps we can follow various
00048 // strategies for this initial 'seed'.
00049 // - start from single cell/face and let FaceCellWave pick up all others
00050 //   from there. might be quite a few waves before everything settles.
00051 // - start from all faces. Lots of initial transfers.
00052 // We do something inbetween:
00053 // - start from all faces where there is a jump. Since we cannot easily
00054 //   determine this across coupled patches (cyclic, processor) introduce
00055 //   all faces of these and let FaceCellWave sort it out.
00056 void smoothDelta::setChangedFaces
00057 (
00058     const polyMesh& mesh,
00059     const volScalarField& delta,
00060     DynamicList<label>& changedFaces,
00061     DynamicList<deltaData>& changedFacesInfo
00062 )
00063 {
00064     for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
00065     {
00066         scalar ownDelta = delta[mesh.faceOwner()[faceI]];
00067 
00068         scalar neiDelta = delta[mesh.faceNeighbour()[faceI]];
00069 
00070         // Check if owner delta much larger than neighbour delta or vice versa
00071 
00072         if (ownDelta > deltaData::maxDeltaRatio * neiDelta)
00073         {
00074             changedFaces.append(faceI);
00075             changedFacesInfo.append(deltaData(ownDelta));
00076         }
00077         else if (neiDelta > deltaData::maxDeltaRatio * ownDelta)
00078         {
00079             changedFaces.append(faceI);
00080             changedFacesInfo.append(deltaData(neiDelta));
00081         }
00082     }
00083 
00084     // Insert all faces of coupled patches no matter what. Let FaceCellWave
00085     // sort it out.
00086     forAll(mesh.boundaryMesh(), patchI) 
00087     {
00088         const polyPatch& patch = mesh.boundaryMesh()[patchI];
00089 
00090         if (patch.coupled())
00091         {
00092             forAll(patch, patchFaceI)
00093             {
00094                 label meshFaceI = patch.start() + patchFaceI;
00095 
00096                 scalar ownDelta = delta[mesh.faceOwner()[meshFaceI]];
00097 
00098                 changedFaces.append(meshFaceI);
00099                 changedFacesInfo.append(deltaData(ownDelta));
00100             }
00101         }
00102     }
00103 
00104     changedFaces.shrink();
00105     changedFacesInfo.shrink();
00106 }
00107 
00108 
00109 void smoothDelta::calcDelta()
00110 {
00111     deltaData::maxDeltaRatio = maxDeltaRatio_;
00112     const volScalarField& geometricDelta = geometricDelta_();
00113 
00114     // Fill changed faces with info
00115     DynamicList<label> changedFaces(mesh_.nFaces()/100 + 100);
00116     DynamicList<deltaData> changedFacesInfo(changedFaces.size());
00117 
00118     setChangedFaces(mesh_, geometricDelta, changedFaces, changedFacesInfo);
00119 
00120     // Set initial field on cells.
00121     List<deltaData> cellDeltaData(mesh_.nCells());
00122 
00123     forAll(geometricDelta, cellI)
00124     {
00125         cellDeltaData[cellI] = geometricDelta[cellI];
00126     }
00127 
00128     // Set initial field on faces.
00129     List<deltaData> faceDeltaData(mesh_.nFaces());
00130 
00131 
00132     // Propagate information over whole domain.
00133     FaceCellWave<deltaData> deltaCalc
00134     (
00135         mesh_,
00136         changedFaces,
00137         changedFacesInfo,
00138         faceDeltaData,
00139         cellDeltaData,
00140         mesh_.globalData().nTotalCells()    // max iterations
00141     );
00142 
00143     forAll(delta_, cellI)
00144     {
00145         delta_[cellI] = cellDeltaData[cellI].delta();
00146     }
00147 }
00148 
00149 
00150 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00151 
00152 smoothDelta::smoothDelta
00153 (
00154     const word& name,
00155     const fvMesh& mesh,
00156     const dictionary& dd
00157 )
00158 :
00159     LESdelta(name, mesh),
00160     geometricDelta_
00161     (
00162         LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
00163     ),
00164     maxDeltaRatio_
00165     (
00166         readScalar(dd.subDict(type() + "Coeffs").lookup("maxDeltaRatio"))
00167     )
00168 {
00169     calcDelta();
00170 }
00171 
00172 
00173 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00174 
00175 void smoothDelta::read(const dictionary& d)
00176 {
00177     const dictionary& dd(d.subDict(type() + "Coeffs"));
00178 
00179     geometricDelta_().read(dd);
00180     dd.lookup("maxDeltaRatio") >> maxDeltaRatio_;
00181     calcDelta();
00182 }
00183 
00184 
00185 void smoothDelta::correct()
00186 {
00187     geometricDelta_().correct();
00188 
00189     if (mesh_.changing())
00190     {
00191         calcDelta();
00192     }
00193 }
00194 
00195 
00196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00197 
00198 } // End namespace Foam
00199 
00200 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines