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

icoMomentError.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     icoMomentError
00026 
00027 Description
00028     Estimates error for the incompressible laminar CFD application icoFoam.
00029 
00030 Usage
00031 
00032     - icoMomentError [OPTIONS]
00033 
00034     @param -noZero \n
00035     Ignore timestep 0.
00036 
00037     @param -constant \n
00038     Include the constant directory.
00039 
00040     @param -time <time>\n
00041     Apply only to specific time.
00042 
00043     @param -latestTime \n
00044     Only apply to latest time step.
00045 
00046     @param -case <dir>\n
00047     Case directory.
00048 
00049     @param -parallel \n
00050     Run in parallel.
00051 
00052     @param -help \n
00053     Display help message.
00054 
00055     @param -doc \n
00056     Display Doxygen API documentation page for this application.
00057 
00058     @param -srcDoc \n
00059     Display Doxygen source documentation page for this application.
00060 
00061 \*---------------------------------------------------------------------------*/
00062 
00063 #include <finiteVolume/fvCFD.H>
00064 #include <finiteVolume/linear.H>
00065 #include <finiteVolume/gaussConvectionScheme.H>
00066 #include <finiteVolume/gaussLaplacianScheme.H>
00067 
00068 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00069 
00070 int main(int argc, char *argv[])
00071 {
00072     timeSelector::addOptions();
00073 
00074 #   include <OpenFOAM/setRootCase.H>
00075 #   include <OpenFOAM/createTime.H>
00076 
00077     instantList timeDirs = timeSelector::select0(runTime, args);
00078 
00079 #   include <OpenFOAM/createMesh.H>
00080 
00081     Info<< "\nEstimating error in the incompressible momentum equation\n"
00082         << "Reading transportProperties\n" << endl;
00083 
00084     IOdictionary transportProperties
00085     (
00086         IOobject
00087         (
00088             "transportProperties",
00089             runTime.constant(),
00090             mesh,
00091             IOobject::MUST_READ,
00092             IOobject::NO_WRITE
00093         )
00094     );
00095 
00096     dimensionedScalar nu
00097     (
00098         transportProperties.lookup("nu")
00099     );
00100 
00101     forAll(timeDirs, timeI)
00102     {
00103         runTime.setTime(timeDirs[timeI], timeI);
00104 
00105         Info<< "Time = " << runTime.timeName() << endl;
00106 
00107         mesh.readUpdate();
00108 
00109         IOobject pHeader
00110         (
00111             "p",
00112             runTime.timeName(),
00113             mesh,
00114             IOobject::MUST_READ
00115         );
00116 
00117         IOobject Uheader
00118         (
00119             "U",
00120             runTime.timeName(),
00121             mesh,
00122             IOobject::MUST_READ
00123         );
00124 
00125         if (pHeader.headerOk() && Uheader.headerOk())
00126         {
00127             Info << "Reading p" << endl;
00128             volScalarField p(pHeader, mesh);
00129 
00130             Info << "Reading U" << endl;
00131             volVectorField U(Uheader, mesh);
00132 
00133 #           include <finiteVolume/createPhi.H>
00134 
00135             volScalarField ek = 0.5*magSqr(U);
00136             volTensorField gradU = fvc::grad(U);
00137 
00138             // Divergence of the error in U squared
00139 
00140             volScalarField L
00141             (
00142                 IOobject
00143                 (
00144                     "L",
00145                     mesh.time().timeName(),
00146                     mesh,
00147                     IOobject::NO_READ,
00148                     IOobject::NO_WRITE
00149                 ),
00150                 mesh,
00151                 dimensionedScalar("one", dimLength, 1.0)
00152             );
00153 
00154             L.internalField() =
00155                 mesh.V()/fvc::surfaceSum(mesh.magSf())().internalField();
00156 
00157             // Warning: 4th row of this equation specially modified
00158             // for the momentum equation. The "real" formulation would
00159             // have diffusivity*(gradV && gradV)
00160             volScalarField momError
00161             (
00162                 IOobject
00163                 (
00164                     "momErrorL" + U.name(),
00165                     mesh.time().timeName(),
00166                     mesh,
00167                     IOobject::NO_READ,
00168                     IOobject::NO_WRITE
00169                 ),
00170                 sqrt
00171                 (
00172                     2.0*mag
00173                     (
00174                         (
00175                             fv::gaussConvectionScheme<scalar>
00176                             (
00177                                 mesh,
00178                                 phi,
00179                                 tmp<surfaceInterpolationScheme<scalar> >
00180                                 (
00181                                     new linear<scalar>(mesh)
00182                                 )
00183                             ).fvcDiv(phi, ek)
00184 
00185                           - nu*
00186                             fv::gaussLaplacianScheme<scalar, scalar>(mesh)
00187                            .fvcLaplacian
00188                             (
00189                                 ek
00190                             )
00191                           - (U & fvc::grad(p))
00192 //                        + nu*(gradU && gradU)
00193                           + 0.5*nu*
00194                             (
00195                                 gradU && (gradU + gradU.T())
00196                             )
00197                         )*L/(mag(U) + nu/L)
00198                     )
00199                 )
00200             );
00201 
00202             momError.boundaryField() = 0.0;
00203             momError.write();
00204         }
00205         else
00206         {
00207             Info<< "    No p or U" << endl;
00208         }
00209 
00210         Info<< endl;
00211     }
00212 
00213     Info<< "End\n" << endl;
00214 
00215     return 0;
00216 }
00217 
00218 
00219 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines