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

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