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

icoErrorEstimate.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     icoErrorEstimate
00026 
00027 Description
00028     Estimates error for the incompressible laminar CFD application icoFoam.
00029 
00030 Usage
00031 
00032     - icoErrorEstimate [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 <errorEstimation/errorEstimate.H>
00065 #include <errorEstimation/resError.H>
00066 
00067 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00068 
00069 int main(int argc, char *argv[])
00070 {
00071     timeSelector::addOptions();
00072 
00073 #   include <OpenFOAM/setRootCase.H>
00074 #   include <OpenFOAM/createTime.H>
00075 
00076     instantList timeDirs = timeSelector::select0(runTime, args);
00077 
00078 #   include <OpenFOAM/createMesh.H>
00079 
00080     Info<< "\nEstimating error in the incompressible momentum equation\n"
00081         << "Reading transportProperties\n" << endl;
00082 
00083     IOdictionary transportProperties
00084     (
00085         IOobject
00086         (
00087             "transportProperties",
00088             runTime.constant(),
00089             mesh,
00090             IOobject::MUST_READ,
00091             IOobject::NO_WRITE
00092         )
00093     );
00094 
00095     dimensionedScalar nu
00096     (
00097         transportProperties.lookup("nu")
00098     );
00099 
00100     forAll(timeDirs, timeI)
00101     {
00102         runTime.setTime(timeDirs[timeI], timeI);
00103 
00104         Info<< "Time = " << runTime.timeName() << endl;
00105 
00106         mesh.readUpdate();
00107 
00108         IOobject pHeader
00109         (
00110             "p",
00111             runTime.timeName(),
00112             mesh,
00113             IOobject::MUST_READ
00114         );
00115 
00116         IOobject Uheader
00117         (
00118             "U",
00119             runTime.timeName(),
00120             mesh,
00121             IOobject::MUST_READ
00122         );
00123 
00124         if (pHeader.headerOk() && Uheader.headerOk())
00125         {
00126             Info << "Reading p" << endl;
00127             volScalarField p(pHeader, mesh);
00128 
00129             Info << "Reading U" << endl;
00130             volVectorField U(Uheader, mesh);
00131 
00132 #           include <finiteVolume/createPhi.H>
00133 
00134             errorEstimate<vector> ee
00135             (
00136                 resError::div(phi, U)
00137               - resError::laplacian(nu, U)
00138              ==
00139                -fvc::grad(p)
00140             );
00141 
00142             volVectorField e = ee.error();
00143             e.write();
00144             mag(e)().write();
00145         }
00146         else
00147         {
00148             Info<< "    No p or U" << endl;
00149         }
00150 
00151         Info<< endl;
00152     }
00153 
00154     Info<< "End\n" << endl;
00155 
00156     return 0;
00157 }
00158 
00159 
00160 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines