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

extendedLeastSquaresVectors.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 "extendedLeastSquaresVectors.H"
00027 #include <finiteVolume/surfaceFields.H>
00028 #include <finiteVolume/volFields.H>
00029 
00030 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00031 
00032 namespace Foam
00033 {
00034     defineTypeNameAndDebug(extendedLeastSquaresVectors, 0);
00035 }
00036 
00037 
00038 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
00039 
00040 Foam::extendedLeastSquaresVectors::extendedLeastSquaresVectors
00041 (
00042     const fvMesh& mesh,
00043     const scalar minDet
00044 )
00045 :
00046     MeshObject<fvMesh, extendedLeastSquaresVectors>(mesh),
00047     minDet_(minDet),
00048     pVectorsPtr_(NULL),
00049     nVectorsPtr_(NULL),
00050     additionalCellsPtr_(NULL),
00051     additionalVectorsPtr_(NULL)
00052 {}
00053 
00054 
00055 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
00056 
00057 Foam::extendedLeastSquaresVectors::~extendedLeastSquaresVectors()
00058 {
00059     deleteDemandDrivenData(pVectorsPtr_);
00060     deleteDemandDrivenData(nVectorsPtr_);
00061 
00062     deleteDemandDrivenData(additionalCellsPtr_);
00063     deleteDemandDrivenData(additionalVectorsPtr_);
00064 }
00065 
00066 
00067 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00068 
00069 void Foam::extendedLeastSquaresVectors::makeLeastSquaresVectors() const
00070 {
00071     if (debug)
00072     {
00073         Info<< "extendedLeastSquaresVectors::makeLeastSquaresVectors() :"
00074             << "Constructing least square gradient vectors"
00075             << endl;
00076     }
00077 
00078     const fvMesh& mesh = mesh_;
00079 
00080     pVectorsPtr_ = new surfaceVectorField
00081     (
00082         IOobject
00083         (
00084             "LeastSquaresP",
00085             mesh_.pointsInstance(),
00086             mesh_,
00087             IOobject::NO_READ,
00088             IOobject::NO_WRITE,
00089             false
00090         ),
00091         mesh_,
00092         dimensionedVector("zero", dimless/dimLength, vector::zero)
00093     );
00094     surfaceVectorField& lsP = *pVectorsPtr_;
00095 
00096     nVectorsPtr_ = new surfaceVectorField
00097     (
00098         IOobject
00099         (
00100             "LeastSquaresN",
00101             mesh_.pointsInstance(),
00102             mesh_,
00103             IOobject::NO_READ,
00104             IOobject::NO_WRITE,
00105             false
00106         ),
00107         mesh_,
00108         dimensionedVector("zero", dimless/dimLength, vector::zero)
00109     );
00110     surfaceVectorField& lsN = *nVectorsPtr_;
00111 
00112     // Set local references to mesh data
00113     const unallocLabelList& owner = mesh_.owner();
00114     const unallocLabelList& neighbour = mesh_.neighbour();
00115 
00116     // Build the d-vectors
00117     surfaceVectorField d = mesh_.Sf()/(mesh_.magSf()*mesh_.deltaCoeffs());
00118 
00119     if (!mesh_.orthogonal())
00120     {
00121         d -= mesh_.correctionVectors()/mesh_.deltaCoeffs();
00122     }
00123 
00124 
00125     // Set up temporary storage for the dd tensor (before inversion)
00126     symmTensorField dd(mesh_.nCells(), symmTensor::zero);
00127 
00128     forAll(owner, faceI)
00129     {
00130         symmTensor wdd = 1.0/magSqr(d[faceI])*sqr(d[faceI]);
00131 
00132         dd[owner[faceI]] += wdd;
00133         dd[neighbour[faceI]] += wdd;
00134     }
00135 
00136     // Visit the boundaries. Coupled boundaries are taken into account
00137     // in the construction of d vectors.  
00138     forAll(d.boundaryField(), patchI)
00139     {
00140         const fvsPatchVectorField& pd = d.boundaryField()[patchI];
00141 
00142         const fvPatch& p = pd.patch();
00143         const unallocLabelList& faceCells = p.faceCells();
00144 
00145         forAll(pd, patchFaceI)
00146         {
00147             dd[faceCells[patchFaceI]] +=
00148                 (1.0/magSqr(pd[patchFaceI]))*sqr(pd[patchFaceI]);
00149         }
00150     }
00151 
00152     scalarField detdd = det(dd);
00153 
00154     Info<< "max(detdd) = " << max(detdd) << endl;
00155     Info<< "min(detdd) = " << min(detdd) << endl;
00156     Info<< "average(detdd) = " << average(detdd) << endl;
00157 
00158     label nAddCells = 0;
00159     label maxNaddCells = 4*detdd.size();
00160     additionalCellsPtr_ = new List<labelPair>(maxNaddCells);
00161     List<labelPair>& additionalCells_ = *additionalCellsPtr_;
00162 
00163     forAll (detdd, i)
00164     {
00165         label count = 0;
00166 
00167         while (++count < 100 && detdd[i] < minDet_)
00168         {
00169             if (nAddCells == maxNaddCells)
00170             {
00171                 FatalErrorIn
00172                 (
00173                     "extendedLeastSquaresVectors::"
00174                     "makeLeastSquaresVectors() const"
00175                 )   << "nAddCells exceeds maxNaddCells"
00176                     << exit(FatalError);
00177             }
00178 
00179             labelList pointLabels = mesh.cells()[i].labels(mesh.faces());
00180 
00181             scalar maxDetddij = 0.0;
00182 
00183             label addCell = -1;
00184 
00185             forAll(pointLabels, j)
00186             {
00187                 forAll(mesh.pointCells()[pointLabels[j]], k)
00188                 {
00189                     label cellj = mesh.pointCells()[pointLabels[j]][k];
00190 
00191                     if (cellj != i)
00192                     {
00193                         if (cellj != -1)
00194                         {
00195                             vector dCij = (mesh.C()[cellj] - mesh.C()[i]);
00196                             
00197                             symmTensor ddij = 
00198                                 dd[i] + (1.0/magSqr(dCij))*sqr(dCij);
00199 
00200                             scalar detddij = det(ddij);
00201 
00202                             if (detddij > maxDetddij)
00203                             {
00204                                 addCell = cellj;
00205                                 maxDetddij = detddij;
00206                             }
00207                         }
00208                     }
00209                 }
00210             }
00211 
00212             if (addCell != -1)
00213             {
00214                 additionalCells_[nAddCells][0] = i;
00215                 additionalCells_[nAddCells++][1] = addCell;
00216                 vector dCij = mesh.C()[addCell] - mesh.C()[i];
00217                 dd[i] += (1.0/magSqr(dCij))*sqr(dCij);
00218                 detdd[i] = det(dd[i]);
00219             }
00220         }
00221     }
00222 
00223     additionalCells_.setSize(nAddCells);
00224 
00225     Info<< "max(detdd) = " << max(detdd) << endl;
00226     Info<< "min(detdd) = " << min(detdd) << endl;
00227     Info<< "average(detdd) = " << average(detdd) << endl;
00228     Info<< "nAddCells/nCells = " << scalar(nAddCells)/mesh.nCells() << endl;
00229 
00230     // Invert the dd tensor
00231     symmTensorField invDd = inv(dd);
00232 
00233 
00234     // Revisit all faces and calculate the lsP and lsN vectors
00235     forAll(owner, faceI)
00236     {
00237         lsP[faceI] =
00238             (1.0/magSqr(d[faceI]))*(invDd[owner[faceI]] & d[faceI]);
00239 
00240         lsN[faceI] =
00241             ((-1.0)/magSqr(d[faceI]))*(invDd[neighbour[faceI]] & d[faceI]);
00242     }
00243 
00244     forAll(lsP.boundaryField(), patchI)
00245     {
00246         const fvsPatchVectorField& pd = d.boundaryField()[patchI];
00247 
00248         fvsPatchVectorField& patchLsP = lsP.boundaryField()[patchI];
00249 
00250         const fvPatch& p = patchLsP.patch();
00251         const unallocLabelList& faceCells = p.faceCells();
00252 
00253         forAll(p, patchFaceI)
00254         {
00255             patchLsP[patchFaceI] =
00256                 (1.0/magSqr(pd[patchFaceI]))
00257                *(invDd[faceCells[patchFaceI]] & pd[patchFaceI]);
00258         }
00259     }
00260 
00261 
00262     additionalVectorsPtr_ = new vectorField(additionalCells_.size());
00263     vectorField& additionalVectors_ = *additionalVectorsPtr_;
00264 
00265     forAll(additionalCells_, i)
00266     {
00267         vector dCij = 
00268             mesh.C()[additionalCells_[i][1]] - mesh.C()[additionalCells_[i][0]];
00269 
00270         additionalVectors_[i] =
00271             (1.0/magSqr(dCij))*(invDd[additionalCells_[i][0]] & dCij);
00272     }
00273 
00274     if (debug)
00275     {
00276         Info<< "extendedLeastSquaresVectors::makeLeastSquaresVectors() :"
00277             << "Finished constructing least square gradient vectors"
00278             << endl;
00279     }
00280 }
00281 
00282 
00283 const Foam::surfaceVectorField&
00284 Foam::extendedLeastSquaresVectors::pVectors() const
00285 {
00286     if (!pVectorsPtr_)
00287     {
00288         makeLeastSquaresVectors();
00289     }
00290 
00291     return *pVectorsPtr_;
00292 }
00293 
00294 
00295 const Foam::surfaceVectorField&
00296 Foam::extendedLeastSquaresVectors::nVectors() const
00297 {
00298     if (!nVectorsPtr_)
00299     {
00300         makeLeastSquaresVectors();
00301     }
00302 
00303     return *nVectorsPtr_;
00304 }
00305 
00306 
00307 const Foam::List<Foam::labelPair>&
00308 Foam::extendedLeastSquaresVectors::additionalCells() const
00309 {
00310     if (!additionalCellsPtr_)
00311     {
00312         makeLeastSquaresVectors();
00313     }
00314 
00315     return *additionalCellsPtr_;
00316 }
00317 
00318 
00319 const Foam::vectorField&
00320 Foam::extendedLeastSquaresVectors::additionalVectors() const
00321 {
00322     if (!additionalVectorsPtr_)
00323     {
00324         makeLeastSquaresVectors();
00325     }
00326 
00327     return *additionalVectorsPtr_;
00328 }
00329 
00330 
00331 bool Foam::extendedLeastSquaresVectors::movePoints()
00332 {
00333     deleteDemandDrivenData(pVectorsPtr_);
00334     deleteDemandDrivenData(nVectorsPtr_);
00335 
00336     deleteDemandDrivenData(additionalCellsPtr_);
00337     deleteDemandDrivenData(additionalVectorsPtr_);
00338 
00339     return true;
00340 }
00341 
00342 
00343 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines