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

smoothDelta.H

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 Class
00025     Foam::smoothDelta
00026 
00027 Description
00028     Smoothed delta which takes a given simple geometric delta and applies
00029     smoothing to it such that the ratio of deltas between two cells is no
00030     larger than a specified amount, typically 1.15.
00031 
00032 SourceFiles
00033     smoothDelta.C
00034 
00035 \*---------------------------------------------------------------------------*/
00036 
00037 #ifndef smoothDelta_H
00038 #define smoothDelta_H
00039 
00040 #include <LESdeltas/LESdelta.H>
00041 
00042 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00043 
00044 namespace Foam
00045 {
00046 
00047 /*---------------------------------------------------------------------------*\
00048                            Class smoothDelta Declaration
00049 \*---------------------------------------------------------------------------*/
00050 
00051 class smoothDelta
00052 :
00053     public LESdelta
00054 {
00055     // Private data
00056 
00057         autoPtr<LESdelta> geometricDelta_;
00058         scalar maxDeltaRatio_;
00059 
00060 
00061     // Private Member Functions
00062 
00063         //- Disallow default bitwise copy construct and assignment
00064         smoothDelta(const smoothDelta&);
00065         void operator=(const smoothDelta&);
00066 
00067         // Calculate the delta values
00068         void calcDelta();
00069 
00070         //- Private member class used by mesh-wave to propagate the delta-ratio
00071         class deltaData
00072         {
00073             scalar delta_;
00074 
00075             // Private Member Functions
00076 
00077                 //- Update. Gets information from neighbouring face/cell and
00078                 //  uses this to update itself (if nessecary) and return true.
00079                 inline bool update
00080                 (
00081                     const deltaData& w2,
00082                     const scalar scale,
00083                     const scalar tol
00084                 );
00085 
00086 
00087         public:
00088 
00089             // Static data members
00090 
00091                 //- delta fraction
00092                 static scalar maxDeltaRatio;
00093 
00094 
00095                 // Constructors
00096 
00097                 //- Construct null
00098                 inline deltaData();
00099 
00100                 //- Construct from origin, yStar, distance
00101                 inline deltaData(const scalar delta);
00102 
00103 
00104             // Member Functions
00105 
00106                 // Access
00107 
00108                 scalar delta() const
00109                 {
00110                     return delta_;
00111                 }
00112 
00113 
00114                 // Needed by FaceCellWave
00115 
00116                     //- Check whether origin has been changed at all or
00117                     //  still contains original (invalid) value.
00118                     inline bool valid() const;
00119 
00120                     //- Check for identical geometrical data.
00121                     //  Used for cyclics checking.
00122                     inline bool sameGeometry
00123                     (
00124                         const polyMesh&,
00125                         const deltaData&,
00126                         const scalar
00127                     ) const;
00128 
00129                     //- Convert any absolute coordinates into relative to
00130                     //  (patch)face centre
00131                     inline void leaveDomain
00132                     (
00133                         const polyMesh&,
00134                         const polyPatch&,
00135                         const label patchFaceI,
00136                         const point& faceCentre
00137                     );
00138 
00139                     //- Reverse of leaveDomain
00140                     inline void enterDomain
00141                     (
00142                         const polyMesh&,
00143                         const polyPatch&,
00144                         const label patchFaceI,
00145                         const point& faceCentre
00146                     );
00147 
00148                     //- Apply rotation matrix to any coordinates
00149                     inline void transform
00150                     (
00151                         const polyMesh&,
00152                         const tensor&
00153                     );
00154 
00155                     //- Influence of neighbouring face.
00156                     inline bool updateCell
00157                     (
00158                         const polyMesh&,
00159                         const label thisCellI,
00160                         const label neighbourFaceI,
00161                         const deltaData& neighbourInfo,
00162                         const scalar tol
00163                     );
00164 
00165                     //- Influence of neighbouring cell.
00166                     inline bool updateFace
00167                     (
00168                         const polyMesh&,
00169                         const label thisFaceI,
00170                         const label neighbourCellI,
00171                         const deltaData& neighbourInfo,
00172                         const scalar tol
00173                     );
00174 
00175                     //- Influence of different value on same face.
00176                     inline bool updateFace
00177                     (
00178                         const polyMesh&,
00179                         const label thisFaceI,
00180                         const deltaData& neighbourInfo,
00181                         const scalar tol
00182                     );
00183 
00184                 // Member Operators
00185 
00186                     // Needed for List IO
00187                     inline bool operator==(const deltaData&) const;
00188 
00189                     inline bool operator!=(const deltaData&) const;
00190 
00191                 // IOstream Operators
00192 
00193                     friend Ostream& operator<<
00194                     (
00195                         Ostream& os,
00196                         const deltaData& wDist
00197                     )
00198                     {
00199                         return os << wDist.delta_;
00200                     }
00201 
00202                     friend Istream& operator>>(Istream& is, deltaData& wDist)
00203                     {
00204                         return is >> wDist.delta_;
00205                     }
00206         };
00207 
00208 
00209         void setChangedFaces
00210         (
00211             const polyMesh& mesh,
00212             const volScalarField& delta,
00213             DynamicList<label>& changedFaces,
00214             DynamicList<deltaData>& changedFacesInfo
00215         );
00216 
00217 
00218 public:
00219 
00220     //- Runtime type information
00221     TypeName("smooth");
00222 
00223 
00224     // Constructors
00225 
00226         //- Construct from name, mesh and IOdictionary
00227         smoothDelta
00228         (
00229             const word& name,
00230             const fvMesh& mesh,
00231             const dictionary&
00232         );
00233 
00234 
00235     //- Destructor
00236     virtual ~smoothDelta()
00237     {}
00238 
00239 
00240     // Member Functions
00241 
00242         //- Read the LESdelta dictionary
00243         virtual void read(const dictionary&);
00244 
00245         // Correct values
00246         virtual void correct();
00247 };
00248 
00249 
00250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00251 
00252 } // End namespace Foam
00253 
00254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00255 
00256 #include <LESdeltas/smoothDeltaDeltaDataI.H>
00257 
00258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00259 
00260 #endif
00261 
00262 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines