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

estimateScalarError.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     estimateScalarError
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     - estimateScalarError [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 <errorEstimation/errorEstimate.H>
00066 #include <errorEstimation/resError.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 scalar transport 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 
00097     Info<< "Reading diffusivity DT\n" << endl;
00098 
00099     dimensionedScalar DT
00100     (
00101         transportProperties.lookup("DT")
00102     );
00103 
00104 
00105     forAll(timeDirs, timeI)
00106     {
00107         runTime.setTime(timeDirs[timeI], timeI);
00108 
00109         Info<< "Time = " << runTime.timeName() << endl;
00110 
00111         mesh.readUpdate();
00112 
00113         IOobject THeader
00114         (
00115             "T",
00116             runTime.timeName(),
00117             mesh,
00118             IOobject::MUST_READ
00119         );
00120 
00121         IOobject Uheader
00122         (
00123             "U",
00124             runTime.timeName(),
00125             mesh,
00126             IOobject::MUST_READ
00127         );
00128 
00129         if (THeader.headerOk() && Uheader.headerOk())
00130         {
00131             Info << "Reading T" << endl;
00132             volScalarField T(THeader, mesh);
00133 
00134             Info << "Reading U" << endl;
00135             volVectorField U(Uheader, mesh);
00136 
00137 #           include <finiteVolume/createPhi.H>
00138 
00139             errorEstimate<scalar> ee
00140             (
00141                 resError::div(phi, T)
00142               - resError::laplacian(DT, T)
00143             );
00144 
00145             ee.residual()().write();
00146             volScalarField e = ee.error();
00147             e.write();
00148             mag(e)().write();
00149         }
00150         else
00151         {
00152             Info<< "    No T or U" << endl;
00153         }
00154 
00155         Info<< endl;
00156     }
00157 
00158     Info<< "End\n" << endl;
00159 
00160     return 0;
00161 }
00162 
00163 
00164 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines